PartMC  2.8.0
run_part.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2017 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_run_part module.
7 
8 !> Monte Carlo simulation.
10 
11  use pmc_constants
12  use pmc_util
13  use pmc_aero_state
14  use pmc_scenario
15  use pmc_env_state
16  use pmc_aero_data
17  use pmc_gas_data
18  use pmc_gas_state
19  use pmc_output
20  use pmc_mosaic
21  use pmc_coagulation
23  use pmc_coag_kernel
24  use pmc_nucleate
25  use pmc_mpi
27  use pmc_photolysis
29 #ifdef PMC_USE_CAMP
30  use camp_camp_core
31  use camp_camp_state
32  use camp_env_state, only: camp_env_state_t => env_state_t
33 #endif
34 #ifdef PMC_USE_SUNDIALS
35  use pmc_condense
36 #endif
37 #ifdef PMC_USE_MPI
38  use mpi
39 #endif
40 
41  !> Type code for undefined or invalid parallel coagulation method.
42  integer, parameter :: parallel_coag_type_invalid = 0
43  !> Type code for local parallel coagulation.
44  integer, parameter :: parallel_coag_type_local = 1
45  !> Type code for distributed parallel coagulation.
46  integer, parameter :: parallel_coag_type_dist = 2
47 
48  !> Options controlling the execution of run_part().
50  !> Final time (s).
51  real(kind=dp) :: t_max
52  !> Output interval (0 disables) (s).
53  real(kind=dp) :: t_output
54  !> Progress interval (0 disables) (s).
55  real(kind=dp) :: t_progress
56  !> Timestep for coagulation.
57  real(kind=dp) :: del_t
58  !> Prefix for output files.
59  character(len=300) :: output_prefix
60  !> Type of coagulation kernel.
61  integer :: coag_kernel_type
62  !> Type of nucleation.
63  integer :: nucleate_type
64  !> Source of nucleation.
65  integer :: nucleate_source
66  !> Weight class of nucleation.
67  integer :: nucleate_weight_class
68  !> Whether to do coagulation.
69  logical :: do_coagulation
70  !> Whether to do nucleation.
71  logical :: do_nucleation
72  !> Allow doubling if needed.
73  logical :: allow_doubling
74  !> Allow halving if needed.
75  logical :: allow_halving
76  !> Whether to do condensation.
77  logical :: do_condensation
78  !> Whether to do MOSAIC.
79  logical :: do_mosaic
80  !> Whether to compute optical properties.
81  logical :: do_optical
82  !> Whether to have explicitly selected weighting.
83  logical :: do_select_weighting
84  !> Type of particle weighting scheme.
85  integer :: weighting_type
86  !> Weighting exponent for power weighting scheme.
87  real(kind=dp) :: weighting_exponent
88  !> Repeat number of run.
89  integer :: i_repeat
90  !> Total number of repeats.
91  integer :: n_repeat
92  !> Wall clock time of start.
93  real(kind=dp) :: t_wall_start
94  !> Whether to record particle removal information.
95  logical :: record_removals
96  !> Whether to run in parallel.
97  logical :: do_parallel
98  !> Parallel output type.
99  integer :: output_type
100  !> Mixing timescale between processes (s).
101  real(kind=dp) :: mix_timescale
102  !> Whether to average gases each timestep.
103  logical :: gas_average
104  !> Whether to average environment each timestep.
105  logical :: env_average
106  !> Parallel coagulation method type.
107  integer :: parallel_coag_type
108  !> Whether to run CAMP.
109  logical :: do_camp_chem
110  !> Whether to run TChem.
111  logical :: do_tchem
112  !> UUID for this simulation.
113  character(len=PMC_UUID_LEN) :: uuid
114  end type run_part_opt_t
115 
116 contains
117 
118 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
119 
120  !> Do a particle-resolved Monte Carlo simulation.
121 #ifdef PMC_USE_CAMP
122  subroutine run_part(scenario, env_state, aero_data, aero_state, gas_data, &
123  gas_state, run_part_opt, camp_core, photolysis)
124 #else
125  subroutine run_part(scenario, env_state, aero_data, aero_state, gas_data, &
126  gas_state, run_part_opt)
127 #endif
128 
129  !> Environment state.
130  type(scenario_t), intent(in) :: scenario
131  !> Environment state.
132  type(env_state_t), intent(inout) :: env_state
133  !> Aerosol data.
134  type(aero_data_t), intent(inout) :: aero_data
135  !> Aerosol state.
136  type(aero_state_t), intent(inout) :: aero_state
137  !> Gas data.
138  type(gas_data_t), intent(in) :: gas_data
139  !> Gas state.
140  type(gas_state_t), intent(inout) :: gas_state
141  !> Monte Carlo options.
142  type(run_part_opt_t), intent(in) :: run_part_opt
143 #ifdef PMC_USE_CAMP
144  !> CAMP chemistry core
145  type(camp_core_t), pointer, intent(inout), optional :: camp_core
146  !> Photolysis calculator
147  type(photolysis_t), pointer, intent(inout), optional :: photolysis
148 #endif
149 
150  real(kind=dp) :: time, pre_time, pre_del_t, prop_done
151  real(kind=dp) :: last_output_time, last_progress_time
152  integer :: rank, n_proc, pre_index, ncid
153  integer :: pre_i_repeat
154  integer :: n_samp, n_coag, n_emit, n_dil_in, n_dil_out, n_nuc
155  integer :: progress_n_samp, progress_n_coag
156  integer :: progress_n_emit, progress_n_dil_in, progress_n_dil_out
157  integer :: progress_n_nuc, n_part_before
158  integer :: global_n_part, global_n_samp, global_n_coag
159  integer :: global_n_emit, global_n_dil_in, global_n_dil_out
160  integer :: global_n_nuc
161  logical :: do_output, do_state, do_state_netcdf, do_progress, did_coag
162  real(kind=dp) :: t_start, t_wall_now, t_wall_elapsed, t_wall_remain
163  type(env_state_t) :: old_env_state
164  integer :: n_time, i_time, pre_i_time
165  integer :: i_state, i_state_netcdf, i_output
166  integer :: i_cur, i_next
167 
168  rank = pmc_mpi_rank()
169  n_proc = pmc_mpi_size()
170 
171  i_time = 0
172  i_output = 1
173  i_state = 1
174  i_state_netcdf = 1
175  time = 0d0
176  progress_n_samp = 0
177  progress_n_coag = 0
178  progress_n_emit = 0
179  progress_n_dil_in = 0
180  progress_n_dil_out = 0
181  progress_n_nuc = 0
182 
183  call check_time_multiple("t_max", run_part_opt%t_max, &
184  "del_t", run_part_opt%del_t)
185  call check_time_multiple("t_output", run_part_opt%t_output, &
186  "del_t", run_part_opt%del_t)
187  call check_time_multiple("t_progress", run_part_opt%t_progress, &
188  "del_t", run_part_opt%del_t)
189 
190  if (run_part_opt%do_mosaic) then
191  call mosaic_init(env_state, aero_data, run_part_opt%del_t, &
192  run_part_opt%do_optical)
193  if (run_part_opt%do_optical) then
194  call mosaic_aero_optical_init(env_state, aero_data, &
195  aero_state, gas_data, gas_state)
196  call mosaic_optical_wavelengths(aero_data)
197  end if
198  end if
199 
200  if (run_part_opt%t_output > 0d0) then
201  call output_state(run_part_opt%output_prefix, &
202  run_part_opt%output_type, aero_data, aero_state, gas_data, &
203  gas_state, env_state, i_state, time, run_part_opt%del_t, &
204  run_part_opt%i_repeat, run_part_opt%record_removals, &
205  run_part_opt%do_optical, run_part_opt%uuid)
206  call aero_info_array_zero(aero_state%aero_info_array)
207  end if
208 
209  call aero_state_rebalance(aero_state, aero_data, &
210  run_part_opt%allow_doubling, &
211  run_part_opt%allow_halving, initial_state_warning=.true.)
212 
213  t_start = env_state%elapsed_time
214  last_output_time = time
215  last_progress_time = time
216  n_time = nint(run_part_opt%t_max / run_part_opt%del_t)
217 
218  global_n_part = aero_state_total_particles_all_procs(aero_state)
219  if (rank == 0) then
220  ! progress only printed from root process
221  if (run_part_opt%i_repeat == 1) then
222  t_wall_elapsed = 0d0
223  t_wall_remain = 0d0
224  else
225  t_wall_now = system_clock_time()
226  prop_done = real(run_part_opt%i_repeat - 1, kind=dp) &
227  / real(run_part_opt%n_repeat, kind=dp)
228  t_wall_elapsed = t_wall_now - run_part_opt%t_wall_start
229  t_wall_remain = (1d0 - prop_done) / prop_done &
230  * t_wall_elapsed
231  end if
232  call print_part_progress(run_part_opt%i_repeat, time, &
233  global_n_part, 0, 0, 0, 0, 0, t_wall_elapsed, t_wall_remain)
234  end if
235 
236  i_cur = 1
237  i_next = n_time
238 
239  call run_part_timeblock(scenario, env_state, aero_data, aero_state, &
240  gas_data, gas_state, run_part_opt, &
241 #ifdef PMC_USE_CAMP
242  camp_core, photolysis, &
243 #endif
244  i_cur, i_next, t_start, last_output_time, last_progress_time, &
245  i_output, progress_n_samp, progress_n_coag, progress_n_emit, &
246  progress_n_dil_in, progress_n_dil_out, progress_n_nuc)
247 
248  if (run_part_opt%do_mosaic) then
249  call mosaic_cleanup()
250  end if
251 
252  end subroutine run_part
253 
254 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
255 
256  !> Print the current simulation progress to the screen.
257  subroutine print_part_progress(i_repeat, t_sim_elapsed, n_part, n_coag, &
258  n_emit, n_dil_in, n_dil_out, n_nuc, t_wall_elapsed, t_wall_remain)
259 
260  !> Repeat number of simulation.
261  integer, intent(in) :: i_repeat
262  !> Elapsed simulation time (s).
263  real(kind=dp), intent(in) :: t_sim_elapsed
264  !> Number of particles.
265  integer, intent(in) :: n_part
266  !> Number of coagulated particles since last progress printing.
267  integer, intent(in) :: n_coag
268  !> Number of emitted particles since last progress printing.
269  integer, intent(in) :: n_emit
270  !> Number of diluted-in particles since last progress printing.
271  integer, intent(in) :: n_dil_in
272  !> Number of diluted-out particles since last progress printing.
273  integer, intent(in) :: n_dil_out
274  !> Number of nucleated particles since last progress printing.
275  integer, intent(in) :: n_nuc
276  !> Elapsed wall time (s).
277  real(kind=dp), intent(in) :: t_wall_elapsed
278  !> Estimated remaining wall time (s).
279  real(kind=dp), intent(in) :: t_wall_remain
280 
281  write(*,'(a6,a1,a6,a1,a7,a1,a7,a1,a7,a1,a8,a1,a9,a1,a7,a1,a6,a1,a6)') &
282  "repeat", " ", "t_sim", " ", "n_part", " ", "n_coag", " ", &
283  "n_emit", " ", "n_dil_in", " ", "n_dil_out", " ", "n_nuc", " ", &
284  "t_wall", " ", "t_rem"
285  write(*,'(a6,a1,a6,a1,a7,a1,a7,a1,a7,a1,a8,a1,a9,a1,a7,a1,a6,a1,a6)') &
286  trim(integer_to_string_max_len(i_repeat, 6)), " ", &
287  trim(time_to_string_max_len(t_sim_elapsed, 6)), " ", &
288  trim(integer_to_string_max_len(n_part, 7)), " ", &
289  trim(integer_to_string_max_len(n_coag, 7)), " ", &
290  trim(integer_to_string_max_len(n_emit, 7)), " ", &
291  trim(integer_to_string_max_len(n_dil_in, 7)), " ", &
292  trim(integer_to_string_max_len(n_dil_out, 7)), " ", &
293  trim(integer_to_string_max_len(n_nuc, 7)), " ", &
294  trim(time_to_string_max_len(t_wall_elapsed, 6)), " ", &
295  trim(time_to_string_max_len(t_wall_remain, 6))
296 
297  end subroutine print_part_progress
298 
299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
300 
301  !> Determines the number of bytes required to pack the given value.
302  integer function pmc_mpi_pack_size_run_part_opt(val)
303 
304  !> Value to pack.
305  type(run_part_opt_t), intent(in) :: val
306 
308  pmc_mpi_pack_size_real(val%t_max) &
309  + pmc_mpi_pack_size_real(val%t_output) &
310  + pmc_mpi_pack_size_real(val%t_progress) &
311  + pmc_mpi_pack_size_real(val%del_t) &
312  + pmc_mpi_pack_size_string(val%output_prefix) &
313  + pmc_mpi_pack_size_integer(val%coag_kernel_type) &
314  + pmc_mpi_pack_size_integer(val%nucleate_type) &
315  + pmc_mpi_pack_size_integer(val%nucleate_source) &
316  + pmc_mpi_pack_size_integer(val%nucleate_weight_class) &
317  + pmc_mpi_pack_size_logical(val%do_coagulation) &
318  + pmc_mpi_pack_size_logical(val%do_nucleation) &
319  + pmc_mpi_pack_size_logical(val%allow_doubling) &
320  + pmc_mpi_pack_size_logical(val%allow_halving) &
321  + pmc_mpi_pack_size_logical(val%do_condensation) &
322  + pmc_mpi_pack_size_logical(val%do_mosaic) &
323  + pmc_mpi_pack_size_logical(val%do_optical) &
324  + pmc_mpi_pack_size_logical(val%do_select_weighting) &
325  + pmc_mpi_pack_size_integer(val%weighting_type) &
326  + pmc_mpi_pack_size_real(val%weighting_exponent) &
327  + pmc_mpi_pack_size_integer(val%i_repeat) &
328  + pmc_mpi_pack_size_integer(val%n_repeat) &
329  + pmc_mpi_pack_size_real(val%t_wall_start) &
330  + pmc_mpi_pack_size_logical(val%record_removals) &
331  + pmc_mpi_pack_size_logical(val%do_parallel) &
332  + pmc_mpi_pack_size_integer(val%output_type) &
333  + pmc_mpi_pack_size_real(val%mix_timescale) &
334  + pmc_mpi_pack_size_logical(val%gas_average) &
335  + pmc_mpi_pack_size_logical(val%env_average) &
336  + pmc_mpi_pack_size_integer(val%parallel_coag_type) &
337  + pmc_mpi_pack_size_logical(val%do_camp_chem) &
338  + pmc_mpi_pack_size_logical(val%do_tchem) &
339  + pmc_mpi_pack_size_string(val%uuid)
340 
341  end function pmc_mpi_pack_size_run_part_opt
342 
343 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
344 
345  !> Packs the given value into the buffer, advancing position.
346  subroutine pmc_mpi_pack_run_part_opt(buffer, position, val)
347 
348  !> Memory buffer.
349  character, intent(inout) :: buffer(:)
350  !> Current buffer position.
351  integer, intent(inout) :: position
352  !> Value to pack.
353  type(run_part_opt_t), intent(in) :: val
354 
355 #ifdef PMC_USE_MPI
356  integer :: prev_position
357 
358  prev_position = position
359  call pmc_mpi_pack_real(buffer, position, val%t_max)
360  call pmc_mpi_pack_real(buffer, position, val%t_output)
361  call pmc_mpi_pack_real(buffer, position, val%t_progress)
362  call pmc_mpi_pack_real(buffer, position, val%del_t)
363  call pmc_mpi_pack_string(buffer, position, val%output_prefix)
364  call pmc_mpi_pack_integer(buffer, position, val%coag_kernel_type)
365  call pmc_mpi_pack_integer(buffer, position, val%nucleate_type)
366  call pmc_mpi_pack_integer(buffer, position, val%nucleate_source)
367  call pmc_mpi_pack_integer(buffer, position, val%nucleate_weight_class)
368  call pmc_mpi_pack_logical(buffer, position, val%do_coagulation)
369  call pmc_mpi_pack_logical(buffer, position, val%do_nucleation)
370  call pmc_mpi_pack_logical(buffer, position, val%allow_doubling)
371  call pmc_mpi_pack_logical(buffer, position, val%allow_halving)
372  call pmc_mpi_pack_logical(buffer, position, val%do_condensation)
373  call pmc_mpi_pack_logical(buffer, position, val%do_mosaic)
374  call pmc_mpi_pack_logical(buffer, position, val%do_optical)
375  call pmc_mpi_pack_logical(buffer, position, val%do_select_weighting)
376  call pmc_mpi_pack_integer(buffer, position, val%weighting_type)
377  call pmc_mpi_pack_real(buffer, position, val%weighting_exponent)
378  call pmc_mpi_pack_integer(buffer, position, val%i_repeat)
379  call pmc_mpi_pack_integer(buffer, position, val%n_repeat)
380  call pmc_mpi_pack_real(buffer, position, val%t_wall_start)
381  call pmc_mpi_pack_logical(buffer, position, val%record_removals)
382  call pmc_mpi_pack_logical(buffer, position, val%do_parallel)
383  call pmc_mpi_pack_integer(buffer, position, val%output_type)
384  call pmc_mpi_pack_real(buffer, position, val%mix_timescale)
385  call pmc_mpi_pack_logical(buffer, position, val%gas_average)
386  call pmc_mpi_pack_logical(buffer, position, val%env_average)
387  call pmc_mpi_pack_integer(buffer, position, val%parallel_coag_type)
388  call pmc_mpi_pack_logical(buffer, position, val%do_camp_chem)
389  call pmc_mpi_pack_logical(buffer, position, val%do_tchem)
390  call pmc_mpi_pack_string(buffer, position, val%uuid)
391  call assert(946070052, &
392  position - prev_position <= pmc_mpi_pack_size_run_part_opt(val))
393 #endif
394 
395  end subroutine pmc_mpi_pack_run_part_opt
396 
397 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
398 
399  !> Unpacks the given value from the buffer, advancing position.
400  subroutine pmc_mpi_unpack_run_part_opt(buffer, position, val)
401 
402  !> Memory buffer.
403  character, intent(inout) :: buffer(:)
404  !> Current buffer position.
405  integer, intent(inout) :: position
406  !> Value to pack.
407  type(run_part_opt_t), intent(inout) :: val
408 
409 #ifdef PMC_USE_MPI
410  integer :: prev_position
411 
412  prev_position = position
413  call pmc_mpi_unpack_real(buffer, position, val%t_max)
414  call pmc_mpi_unpack_real(buffer, position, val%t_output)
415  call pmc_mpi_unpack_real(buffer, position, val%t_progress)
416  call pmc_mpi_unpack_real(buffer, position, val%del_t)
417  call pmc_mpi_unpack_string(buffer, position, val%output_prefix)
418  call pmc_mpi_unpack_integer(buffer, position, val%coag_kernel_type)
419  call pmc_mpi_unpack_integer(buffer, position, val%nucleate_type)
420  call pmc_mpi_unpack_integer(buffer, position, val%nucleate_source)
421  call pmc_mpi_unpack_integer(buffer, position, val%nucleate_weight_class)
422  call pmc_mpi_unpack_logical(buffer, position, val%do_coagulation)
423  call pmc_mpi_unpack_logical(buffer, position, val%do_nucleation)
424  call pmc_mpi_unpack_logical(buffer, position, val%allow_doubling)
425  call pmc_mpi_unpack_logical(buffer, position, val%allow_halving)
426  call pmc_mpi_unpack_logical(buffer, position, val%do_condensation)
427  call pmc_mpi_unpack_logical(buffer, position, val%do_mosaic)
428  call pmc_mpi_unpack_logical(buffer, position, val%do_optical)
429  call pmc_mpi_unpack_logical(buffer, position, val%do_select_weighting)
430  call pmc_mpi_unpack_integer(buffer, position, val%weighting_type)
431  call pmc_mpi_unpack_real(buffer, position, val%weighting_exponent)
432  call pmc_mpi_unpack_integer(buffer, position, val%i_repeat)
433  call pmc_mpi_unpack_integer(buffer, position, val%n_repeat)
434  call pmc_mpi_unpack_real(buffer, position, val%t_wall_start)
435  call pmc_mpi_unpack_logical(buffer, position, val%record_removals)
436  call pmc_mpi_unpack_logical(buffer, position, val%do_parallel)
437  call pmc_mpi_unpack_integer(buffer, position, val%output_type)
438  call pmc_mpi_unpack_real(buffer, position, val%mix_timescale)
439  call pmc_mpi_unpack_logical(buffer, position, val%gas_average)
440  call pmc_mpi_unpack_logical(buffer, position, val%env_average)
441  call pmc_mpi_unpack_integer(buffer, position, val%parallel_coag_type)
442  call pmc_mpi_unpack_logical(buffer, position, val%do_camp_chem)
443  call pmc_mpi_unpack_logical(buffer, position, val%do_tchem)
444  call pmc_mpi_unpack_string(buffer, position, val%uuid)
445  call assert(480118362, &
446  position - prev_position <= pmc_mpi_pack_size_run_part_opt(val))
447 #endif
448 
449  end subroutine pmc_mpi_unpack_run_part_opt
450 
451 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
452 
453  !> Read the specification for a run_part simulation from a spec file.
454  subroutine spec_file_read_run_part(file, run_part_opt, aero_data, &
455  aero_state_init, gas_data, gas_state_init, env_state_init, &
456  aero_dist_init, scenario, &
457 #ifdef PMC_USE_CAMP
458  camp_core, photolysis, aero_state, &
459 #endif
460  n_part, rand_init, do_init_equilibrate, do_restart)
461 
462  !> Spec file.
463  type(spec_file_t), intent(inout) :: file
464  !> Monte Carlo options.
465  type(run_part_opt_t), intent(inout) :: run_part_opt
466  !> Aerosol data.
467  type(aero_data_t), intent(inout) :: aero_data
468  !> Initial aerosol state.
469  type(aero_state_t), intent(inout) :: aero_state_init
470  !> Gas data.
471  type(gas_data_t), intent(inout) :: gas_data
472  !> Initial gas state.
473  type(gas_state_t), intent(inout) :: gas_state_init
474  !> Initial environmental state.
475  type(env_state_t), intent(inout) :: env_state_init
476  !> Initial aerosol distribution.
477  type(aero_dist_t), intent(inout) :: aero_dist_init
478  !> Scenario data.
479  type(scenario_t), intent(inout) :: scenario
480 #ifdef PMC_USE_CAMP
481  !> CAMP core.
482  type(camp_core_t), pointer :: camp_core
483  !> Photolysis calculator.
484  type(photolysis_t), pointer :: photolysis
485  !> Aerosol state.
486  type(aero_state_t), intent(inout) :: aero_state
487 #endif
488  !> Ideal number of computational particles.
489  real(kind=dp), intent(inout) :: n_part
490  !> Random number generator seed.
491  integer, intent(out) :: rand_init
492  !> Whether to equilibrate.
493  logical, intent(out) :: do_init_equilibrate
494  !> Whether simulation is a restart.
495  logical, intent(out) :: do_restart
496 
497  integer :: i_repeat, i_group
498  logical :: read_aero_weight_classes
499  character(len=PMC_MAX_FILENAME_LEN) :: restart_filename
500  integer :: dummy_index, dummy_i_repeat
501  real(kind=dp) :: dummy_time, dummy_del_t
502  character(len=PMC_MAX_FILENAME_LEN) :: sub_filename
503  type(spec_file_t) :: sub_file
504  character(len=PMC_MAX_FILENAME_LEN) :: camp_config_filename, &
505  tchem_gas_filename, tchem_aero_filename, tchem_numerics_filename
506  integer :: n_grid_cells
507 
508  call spec_file_read_string(file, 'output_prefix', &
509  run_part_opt%output_prefix)
510  call spec_file_read_integer(file, 'n_repeat', run_part_opt%n_repeat)
511  call spec_file_read_real(file, 'n_part', n_part)
512  call spec_file_read_logical(file, 'restart', do_restart)
513  if (do_restart) then
514  call spec_file_read_string(file, 'restart_file', restart_filename)
515  end if
516 
517  if (.not. do_restart) then
518  call spec_file_read_logical(file, 'do_select_weighting', &
519  run_part_opt%do_select_weighting)
520  read_aero_weight_classes = .false.
521  if (run_part_opt%do_select_weighting) then
522  call spec_file_read_aero_state_weighting_type(file, &
523  run_part_opt%weighting_type, run_part_opt%weighting_exponent)
524  if (run_part_opt%weighting_type &
526  read_aero_weight_classes = .true.
527  end if
528  else
529  run_part_opt%weighting_type = aero_state_weight_nummass_source
530  run_part_opt%weighting_exponent = 0.0d0
531  end if
532  end if
533 
534  call spec_file_read_real(file, 't_max', run_part_opt%t_max)
535  call spec_file_read_real(file, 'del_t', run_part_opt%del_t)
536  call spec_file_read_real(file, 't_output', run_part_opt%t_output)
537  call spec_file_read_real(file, 't_progress', run_part_opt%t_progress)
538 
539  call spec_file_read_logical(file, 'do_camp_chem', run_part_opt%do_camp_chem)
540  if (run_part_opt%do_camp_chem) then
541 #ifdef PMC_USE_CAMP
542  call spec_file_read_string(file, 'camp_config', camp_config_filename)
543  camp_core => camp_core_t(camp_config_filename)
544  call camp_core%initialize()
545  photolysis => photolysis_t(camp_core)
546 #else
547  call spec_file_die_msg(648994111, file, &
548  'cannot do camp chem, CAMP support not compiled in')
549 #endif
550  end if
551 
552  call spec_file_read_logical(file, 'do_tchem', run_part_opt%do_tchem)
553  if (run_part_opt%do_tchem) then
554 #ifdef PMC_USE_TCHEM
555  call spec_file_read_string(file, 'tchem_gas_config', &
556  tchem_gas_filename)
557  call spec_file_read_string(file, 'tchem_aero_config', &
558  tchem_aero_filename)
559  call spec_file_read_string(file, 'tchem_numerics_config', &
560  tchem_numerics_filename)
561 #endif
562  end if
563 
564  if (run_part_opt%do_tchem) then
565 #ifdef PMC_USE_TCHEM
566  n_grid_cells = 1
567  call pmc_tchem_initialize(tchem_gas_filename, tchem_aero_filename, &
568  tchem_numerics_filename, gas_data, aero_data, n_grid_cells)
569 #endif
570  end if
571 
572  if (do_restart) then
573  call input_state(restart_filename, dummy_index, dummy_time, &
574  dummy_del_t, dummy_i_repeat, run_part_opt%uuid, aero_data, &
575  aero_state_init, gas_data, gas_state_init, env_state_init)
576  end if
577 
578  if (.not. do_restart) then
579  env_state_init%elapsed_time = 0d0
580 
581  if (.not. (run_part_opt%do_camp_chem .or. run_part_opt%do_tchem)) then
582  call spec_file_read_string(file, 'gas_data', sub_filename)
583  call spec_file_open(sub_filename, sub_file)
584  call spec_file_read_gas_data(sub_file, gas_data)
585  call spec_file_close(sub_file)
586  else if (run_part_opt%do_camp_chem) then
587 #ifdef PMC_USE_CAMP
588  call gas_data_initialize(gas_data, camp_core)
589 #endif
590  end if
591  call spec_file_read_string(file, 'gas_init', sub_filename)
592  call spec_file_open(sub_filename, sub_file)
593  call spec_file_read_gas_state(sub_file, gas_data, gas_state_init)
594  call spec_file_close(sub_file)
595 
596  if (.not. run_part_opt%do_camp_chem) then
597  call spec_file_read_string(file, 'aerosol_data', sub_filename)
598  call spec_file_open(sub_filename, sub_file)
599  call spec_file_read_aero_data(sub_file, aero_data)
600  call spec_file_close(sub_file)
601  ! FIXME: Temporary to run PartMC. Replace with initialization from TChem
602  else if (run_part_opt%do_tchem) then
603  call spec_file_read_string(file, 'aerosol_data', sub_filename)
604  call spec_file_open(sub_filename, sub_file)
605  call spec_file_read_aero_data(sub_file, aero_data)
606  call spec_file_close(sub_file)
607  else
608 #ifdef PMC_USE_CAMP
609  call aero_data_initialize(aero_data, camp_core)
610  call aero_state_initialize(aero_state, aero_data, camp_core)
611 #endif
612  end if
613  call spec_file_read_fractal(file, aero_data%fractal)
614 
615  call spec_file_read_string(file, 'aerosol_init', sub_filename)
616  call spec_file_open(sub_filename, sub_file)
617  call spec_file_read_aero_dist(sub_file, aero_data, &
618  read_aero_weight_classes, aero_dist_init)
619  call spec_file_close(sub_file)
620  end if
621 
622  call spec_file_read_scenario(file, gas_data, aero_data, &
623  read_aero_weight_classes, scenario)
624  call spec_file_read_env_state(file, env_state_init)
625 
626  call spec_file_read_logical(file, 'do_coagulation', &
627  run_part_opt%do_coagulation)
628  if (run_part_opt%do_coagulation) then
629  call spec_file_read_coag_kernel_type(file, run_part_opt%coag_kernel_type)
630  if (run_part_opt%coag_kernel_type == coag_kernel_type_additive) then
631  call spec_file_read_real(file, 'additive_kernel_coeff', &
632  env_state_init%additive_kernel_coefficient)
633  end if
634  else
635  run_part_opt%coag_kernel_type = coag_kernel_type_invalid
636  end if
637 
638  call spec_file_read_logical(file, 'do_condensation', &
639  run_part_opt%do_condensation)
640 #ifndef PMC_USE_SUNDIALS
641  call assert_msg(121370218, &
642  run_part_opt%do_condensation .eqv. .false., &
643  "cannot use condensation, SUNDIALS support is not compiled in")
644 #endif
645  do_init_equilibrate = .false.
646  if (run_part_opt%do_condensation) then
647  call spec_file_read_logical(file, 'do_init_equilibrate', &
648  do_init_equilibrate)
649  end if
650 
651  call spec_file_read_logical(file, 'do_mosaic', run_part_opt%do_mosaic)
652  if (run_part_opt%do_mosaic .and. (.not. mosaic_support())) then
653  call spec_file_die_msg(230495365, file, &
654  'cannot use MOSAIC, support is not compiled in')
655  end if
656  if (run_part_opt%do_mosaic .and. run_part_opt%do_condensation) then
657  call spec_file_die_msg(599877804, file, &
658  'cannot use MOSAIC and condensation simultaneously')
659  end if
660  if (run_part_opt%do_mosaic) then
661  call spec_file_read_logical(file, 'do_optical', &
662  run_part_opt%do_optical)
663  else
664  run_part_opt%do_optical = .false.
665  end if
666 
667  call spec_file_read_logical(file, 'do_nucleation', &
668  run_part_opt%do_nucleation)
669  if (run_part_opt%do_nucleation) then
670  call spec_file_read_nucleate_type(file, aero_data, &
671  run_part_opt%nucleate_type, run_part_opt%nucleate_source, &
672  run_part_opt%nucleate_weight_class)
673  else
674  run_part_opt%nucleate_type = nucleate_type_invalid
675  end if
676 
677  call spec_file_read_integer(file, 'rand_init', rand_init)
678  call spec_file_read_logical(file, 'allow_doubling', &
679  run_part_opt%allow_doubling)
680  call spec_file_read_logical(file, 'allow_halving', &
681  run_part_opt%allow_halving)
682  call spec_file_read_logical(file, 'record_removals', &
683  run_part_opt%record_removals)
684 
685  call spec_file_read_logical(file, 'do_parallel', run_part_opt%do_parallel)
686  if (run_part_opt%do_parallel) then
687 #ifndef PMC_USE_MPI
688  call spec_file_die_msg(929006383, file, &
689  'cannot use parallel mode, support is not compiled in')
690 #endif
691  call spec_file_read_output_type(file, run_part_opt%output_type)
692  call spec_file_read_real(file, 'mix_timescale', &
693  run_part_opt%mix_timescale)
694  call spec_file_read_logical(file, 'gas_average', &
695  run_part_opt%gas_average)
696  call spec_file_read_logical(file, 'env_average', &
697  run_part_opt%env_average)
698  call spec_file_read_parallel_coag_type(file, &
699  run_part_opt%parallel_coag_type)
700  else
701  run_part_opt%output_type = output_type_single
702  run_part_opt%mix_timescale = 0d0
703  run_part_opt%gas_average = .false.
704  run_part_opt%env_average = .false.
705  run_part_opt%parallel_coag_type = parallel_coag_type_local
706  end if
707 
708  call spec_file_close(file)
709 
710  end subroutine spec_file_read_run_part
711 
712 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
713 
714  !> Read the specification for a parallel coagulation type from a spec file.
715  subroutine spec_file_read_parallel_coag_type(file, parallel_coag_type)
716 
717  !> Spec file.
718  type(spec_file_t), intent(inout) :: file
719  !> Kernel type.
720  integer, intent(out) :: parallel_coag_type
721 
722  character(len=SPEC_LINE_MAX_VAR_LEN) :: parallel_coag_type_name
723 
724  !> \page input_format_parallel_coag Input File Format: Parallel Coagulation Type
725  !!
726  !! The output type is specified by the parameter:
727  !! - \b parallel_coag (string): type of parallel coagulation ---
728  !! must be one of: \c local for only within-process
729  !! coagulation or \c dist to have all processes perform
730  !! coagulation globally, requesting particles from other
731  !! processes as needed
732  !!
733  !! See also:
734  !! - \ref spec_file_format --- the input file text format
735 
736  call spec_file_read_string(file, 'parallel_coag', &
737  parallel_coag_type_name)
738  if (trim(parallel_coag_type_name) == 'local') then
739  parallel_coag_type = parallel_coag_type_local
740  elseif (trim(parallel_coag_type_name) == 'dist') then
741  parallel_coag_type = parallel_coag_type_dist
742  else
743  call spec_file_die_msg(494684716, file, &
744  "Unknown parallel coagulation type: " &
745  // trim(parallel_coag_type_name))
746  end if
747 
748  end subroutine spec_file_read_parallel_coag_type
749 
750 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
751 
752  !> Do a one particle-resolved Monte Carlo simulation timestep.
753  subroutine run_part_timestep(scenario, env_state, aero_data, aero_state, &
754  gas_data, gas_state, run_part_opt, &
755 #ifdef PMC_USE_CAMP
756  camp_core, photolysis, &
757 #endif
758  i_time, t_start, last_output_time, last_progress_time, i_output, &
759  progress_n_samp, progress_n_coag, progress_n_emit, progress_n_dil_in, &
760  progress_n_dil_out, progress_n_nuc)
761 
762  !> Environment state.
763  type(scenario_t), intent(in) :: scenario
764  !> Environment state.
765  type(env_state_t), intent(inout) :: env_state
766  !> Aerosol data.
767  type(aero_data_t), intent(in) :: aero_data
768  !> Aerosol state.
769  type(aero_state_t), intent(inout) :: aero_state
770  !> Gas data.
771  type(gas_data_t), intent(in) :: gas_data
772  !> Gas state.
773  type(gas_state_t), intent(inout) :: gas_state
774  !> Monte Carlo options.
775  type(run_part_opt_t), intent(in) :: run_part_opt
776 #ifdef PMC_USE_CAMP
777  !> CAMP chemistry core
778  type(camp_core_t), pointer, intent(inout), optional :: camp_core
779  !> Photolysis calculator
780  type(photolysis_t), pointer, intent(inout), optional :: photolysis
781 #endif
782  !> Current simulation time.
783  integer, intent(in) :: i_time
784  ! Start time of simulation.
785  real(kind=dp), intent(in) :: t_start
786  !> Last time output was written (s).
787  real(kind=dp), intent(inout) :: last_output_time
788  !> Last time progress was output to screen (s).
789  real(kind=dp), intent(inout) :: last_progress_time
790  !> Output timestep integer for output filename.
791  integer, intent(inout) :: i_output
792  !> Total number of samples tested.
793  integer, intent(inout) :: progress_n_samp
794  !> Number of coagulation events.
795  integer, intent(inout) :: progress_n_coag
796  !> Number of particles added by emission.
797  integer, intent(inout) :: progress_n_emit
798  !> Number of particles added by dilution.
799  integer, intent(inout) :: progress_n_dil_in
800  !> Number of particles removed by dilution.
801  integer, intent(inout) :: progress_n_dil_out
802  !> Number of particles added by nucleation.
803  integer, intent(inout) :: progress_n_nuc
804 
805  real(kind=dp) :: prop_done
806  integer :: n_samp, n_coag, n_emit, n_dil_in, n_dil_out, n_nuc
807  integer :: n_part_before
808  integer :: global_n_part, global_n_samp, global_n_coag
809  integer :: global_n_emit, global_n_dil_in, global_n_dil_out
810  integer :: global_n_nuc
811  logical :: do_output, do_state, do_state_netcdf, do_progress
812  real(kind=dp) :: t_wall_now, t_wall_elapsed, t_wall_remain, time
813  type(env_state_t) :: old_env_state
814 #ifdef PMC_USE_CAMP
815  type(camp_state_t), pointer :: camp_state
816  type(camp_state_t), pointer :: camp_pre_aero_state, camp_post_aero_state
817  type(camp_env_state_t) :: camp_env_state
818 #endif
819 
820  time = i_time * run_part_opt%del_t
821 
822 #ifdef PMC_USE_CAMP
823  if (run_part_opt%do_camp_chem) then
824  camp_env_state%temp = env_state%temp
825  camp_env_state%pressure = env_state%pressure
826  camp_state => camp_core%new_state_one_cell(camp_env_state)
827  camp_pre_aero_state => camp_core%new_state_one_cell(camp_env_state)
828  camp_post_aero_state => camp_core%new_state_one_cell(camp_env_state)
829  end if
830 #endif
831 
832  old_env_state = env_state
833  call scenario_update_env_state(scenario, env_state, time + t_start)
834 
835  if (run_part_opt%do_nucleation) then
836  n_part_before = aero_state_total_particles(aero_state)
837  call nucleate(run_part_opt%nucleate_type, &
838  run_part_opt%nucleate_source, run_part_opt%nucleate_weight_class, &
839  env_state, gas_data, aero_data, &
840  aero_state, gas_state, run_part_opt%del_t, &
841  run_part_opt%allow_doubling, run_part_opt%allow_halving)
842  n_nuc = aero_state_total_particles(aero_state) &
843  - n_part_before
844  progress_n_nuc = progress_n_nuc + n_nuc
845  end if
846 
847  if (run_part_opt%do_coagulation) then
848  if (run_part_opt%parallel_coag_type &
849  == parallel_coag_type_local) then
850  call mc_coag(run_part_opt%coag_kernel_type, env_state, &
851  aero_data, aero_state, run_part_opt%del_t, n_samp, n_coag)
852  elseif (run_part_opt%parallel_coag_type &
853  == parallel_coag_type_dist) then
854  call mc_coag_dist(run_part_opt%coag_kernel_type, env_state, &
855  aero_data, aero_state, run_part_opt%del_t, n_samp, n_coag)
856  else
857  call die_msg(323011762, "unknown parallel coagulation type: " &
858  // trim(integer_to_string(run_part_opt%parallel_coag_type)))
859  end if
860  progress_n_samp = progress_n_samp + n_samp
861  progress_n_coag = progress_n_coag + n_coag
862  end if
863 
864 #ifdef PMC_USE_SUNDIALS
865  if (run_part_opt%do_condensation) then
866  call condense_particles(aero_state, aero_data, old_env_state, &
867  env_state, run_part_opt%del_t)
868  end if
869 #endif
870 
871  call scenario_update_gas_state(scenario, run_part_opt%del_t, &
872  env_state, old_env_state, gas_data, gas_state)
873  call scenario_update_aero_state(scenario, run_part_opt%del_t, &
874  env_state, old_env_state, aero_data, aero_state, n_emit, &
875  n_dil_in, n_dil_out, run_part_opt%allow_doubling, &
876  run_part_opt%allow_halving)
877  progress_n_emit = progress_n_emit + n_emit
878  progress_n_dil_in = progress_n_dil_in + n_dil_in
879  progress_n_dil_out = progress_n_dil_out + n_dil_out
880 
881  if (run_part_opt%do_camp_chem) then
882 #ifdef PMC_USE_CAMP
883  call pmc_camp_interface_solve(camp_core, camp_state, &
884  camp_pre_aero_state, camp_post_aero_state, env_state, &
885  aero_data, aero_state, gas_data, gas_state, photolysis, &
886  run_part_opt%del_t)
887 #endif
888  end if
889 
890  if (run_part_opt%do_tchem) then
891 #ifdef PMC_USE_TCHEM
892  call pmc_tchem_interface_solve(env_state, aero_data, aero_state, &
893  gas_data, gas_state, run_part_opt%del_t)
894 #endif
895  end if
896 
897  if (run_part_opt%do_mosaic) then
898  call mosaic_timestep(env_state, aero_data, aero_state, gas_data, &
899  gas_state, run_part_opt%do_optical, run_part_opt%uuid)
900  end if
901 
902  if (run_part_opt%mix_timescale > 0d0) then
903  call aero_state_mix(aero_state, run_part_opt%del_t, &
904  run_part_opt%mix_timescale, aero_data)
905  end if
906  if (run_part_opt%gas_average) then
907  call gas_state_mix(gas_state)
908  end if
909  if (run_part_opt%gas_average) then
910  call env_state_mix(env_state)
911  end if
912 
913  call aero_state_rebalance(aero_state, aero_data, &
914  run_part_opt%allow_doubling, &
915  run_part_opt%allow_halving, initial_state_warning=.false.)
916 
917 
918  if (run_part_opt%t_output > 0d0) then
919  call check_event(time, run_part_opt%del_t, run_part_opt%t_output, &
920  last_output_time, do_output)
921  if (do_output) then
922  i_output = i_output + 1
923  call output_state(run_part_opt%output_prefix, &
924  run_part_opt%output_type, aero_data, aero_state, gas_data, &
925  gas_state, env_state, i_output, time, run_part_opt%del_t, &
926  run_part_opt%i_repeat, run_part_opt%record_removals, &
927  run_part_opt%do_optical, run_part_opt%uuid)
928  call aero_info_array_zero(aero_state%aero_info_array)
929  end if
930  end if
931 
932  if (.not. run_part_opt%record_removals) then
933  ! If we are not recording removals then we can zero them as
934  ! often as possible to minimize the cost of maintaining
935  ! them.
936  call aero_info_array_zero(aero_state%aero_info_array)
937  end if
938 
939  if (run_part_opt%t_progress > 0d0) then
940  call check_event(time, run_part_opt%del_t, &
941  run_part_opt%t_progress, last_progress_time, do_progress)
942  if (do_progress) then
943  global_n_part = aero_state_total_particles_all_procs(aero_state)
944  call pmc_mpi_reduce_sum_integer(progress_n_samp, global_n_samp)
945  call pmc_mpi_reduce_sum_integer(progress_n_coag, global_n_coag)
946  call pmc_mpi_reduce_sum_integer(progress_n_emit, global_n_emit)
947  call pmc_mpi_reduce_sum_integer(progress_n_dil_in, &
948  global_n_dil_in)
949  call pmc_mpi_reduce_sum_integer(progress_n_dil_out, &
950  global_n_dil_out)
951  call pmc_mpi_reduce_sum_integer(progress_n_nuc, global_n_nuc)
952  if (pmc_mpi_rank() == 0) then
953  ! progress only printed from root process
954  t_wall_now = system_clock_time()
955  prop_done = (real(run_part_opt%i_repeat - 1, kind=dp) &
956  + time / run_part_opt%t_max) &
957  / real(run_part_opt%n_repeat, kind=dp)
958  t_wall_elapsed = t_wall_now - run_part_opt%t_wall_start
959  t_wall_remain = (1d0 - prop_done) / prop_done * t_wall_elapsed
960  call print_part_progress(run_part_opt%i_repeat, time, &
961  global_n_part, global_n_coag, global_n_emit, &
962  global_n_dil_in, global_n_dil_out, global_n_nuc, &
963  t_wall_elapsed, t_wall_remain)
964  end if
965  ! reset counters so they show information since last
966  ! progress display
967  progress_n_samp = 0
968  progress_n_coag = 0
969  progress_n_emit = 0
970  progress_n_dil_in = 0
971  progress_n_dil_out = 0
972  progress_n_nuc = 0
973  end if
974  end if
975 
976  end subroutine run_part_timestep
977 
978 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
979 
980  !> Do a number of time steps of particle-reoslved Monte Carlo simulation.
981  subroutine run_part_timeblock(scenario, env_state, aero_data, aero_state, &
982  gas_data, gas_state, run_part_opt, &
983 #ifdef PMC_USE_CAMP
984  camp_core, photolysis, &
985 #endif
986  i_cur, i_next, t_start, last_output_time, last_progress_time, &
987  i_output, progress_n_samp, progress_n_coag, progress_n_emit, &
988  progress_n_dil_in, progress_n_dil_out, progress_n_nuc)
989 
990  !> Environment state.
991  type(scenario_t), intent(in) :: scenario
992  !> Environment state.
993  type(env_state_t), intent(inout) :: env_state
994  !> Aerosol data.
995  type(aero_data_t), intent(in) :: aero_data
996  !> Aerosol state.
997  type(aero_state_t), intent(inout) :: aero_state
998  !> Gas data.
999  type(gas_data_t), intent(in) :: gas_data
1000  !> Gas state.
1001  type(gas_state_t), intent(inout) :: gas_state
1002  !> Monte Carlo options.
1003  type(run_part_opt_t), intent(in) :: run_part_opt
1004 #ifdef PMC_USE_CAMP
1005  !> CAMP chemistry core
1006  type(camp_core_t), pointer, intent(inout), optional :: camp_core
1007  !> Photolysis calculator
1008  type(photolysis_t), pointer, intent(inout), optional :: photolysis
1009 #endif
1010  !> Current simulation timestep.
1011  integer, intent(in) :: i_cur
1012  ! Start time of simulation.
1013  real(kind=dp), intent(in) :: t_start
1014  ! End timestep to simulate.
1015  integer, intent(in) :: i_next
1016  !> Last time output was written (s).
1017  real(kind=dp), intent(inout) :: last_output_time
1018  !> Last time progress was output to screen (s).
1019  real(kind=dp), intent(inout) :: last_progress_time
1020  !> Output timestep integer for output filename.
1021  integer, intent(inout) :: i_output
1022  !> Total number of samples tested.
1023  integer, intent(inout) :: progress_n_samp
1024  !> Number of coagulation events.
1025  integer, intent(inout) :: progress_n_coag
1026  !> Number of particles added by emission.
1027  integer, intent(inout) :: progress_n_emit
1028  !> Number of particles added by dilution.
1029  integer, intent(inout) :: progress_n_dil_in
1030  !> Number of particles removed by dilution.
1031  integer, intent(inout) :: progress_n_dil_out
1032  !> Number of particles added by nucleation.
1033  integer, intent(inout) :: progress_n_nuc
1034 
1035  integer :: i_time
1036 
1037  do i_time = i_cur,i_next
1038  call run_part_timestep(scenario, env_state, aero_data, aero_state, &
1039  gas_data, gas_state, run_part_opt, &
1040 #ifdef PMC_USE_CAMP
1041  camp_core, photolysis, &
1042 #endif
1043  i_time, t_start, last_output_time, last_progress_time, i_output, &
1044  progress_n_samp, progress_n_coag, progress_n_emit, &
1045  progress_n_dil_in, progress_n_dil_out, progress_n_nuc)
1046  end do
1047 
1048  end subroutine run_part_timeblock
1049 
1050 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1051 
1052  !> Read the specification for a run_part simulation from a spec file.
1053  subroutine pmc_mpi_broadcast_run_part(run_part_opt, aero_data, &
1054  aero_state_init, gas_data, gas_state_init, env_state_init, &
1055  aero_dist_init, scenario, &
1056 #ifdef PMC_USE_CAMP
1057  camp_core, photolysis, aero_state, &
1058 #endif
1059  n_part, rand_init, do_init_equilibrate, do_restart)
1060 
1061  !> Monte Carlo options.
1062  type(run_part_opt_t), intent(inout) :: run_part_opt
1063  !> Aerosol data.
1064  type(aero_data_t), intent(inout) :: aero_data
1065  !> Initial aerosol state.
1066  type(aero_state_t), intent(inout) :: aero_state_init
1067  !> Gas data.
1068  type(gas_data_t), intent(inout) :: gas_data
1069  !> Initial gas state.
1070  type(gas_state_t), intent(inout) :: gas_state_init
1071  !> Initial environmental state.
1072  type(env_state_t), intent(inout) :: env_state_init
1073  !> Initial aerosol distribution.
1074  type(aero_dist_t), intent(inout) :: aero_dist_init
1075  !> Scenario data.
1076  type(scenario_t), intent(inout) :: scenario
1077 #ifdef PMC_USE_CAMP
1078  !> CAMP core.
1079  type(camp_core_t), pointer :: camp_core
1080  !> Photolysis calculator.
1081  type(photolysis_t), pointer :: photolysis
1082  !> Aerosol state.
1083  type(aero_state_t), intent(inout) :: aero_state
1084 #endif
1085  !> Ideal number of computational particles.
1086  real(kind=dp), intent(inout) :: n_part
1087  !> Random number generator seed.
1088  integer, intent(inout) :: rand_init
1089  !> Whether to equilibrate.
1090  logical, intent(inout) :: do_init_equilibrate
1091  !> Whether simulation is a restart.
1092  logical, intent(inout) :: do_restart
1093 
1094  character, allocatable :: buffer(:)
1095  integer :: buffer_size, max_buffer_size
1096  integer :: position
1097 
1098 #ifdef PMC_USE_MPI
1099  if (pmc_mpi_rank() == 0) then
1100  ! root process determines size
1101  max_buffer_size = 0
1102  max_buffer_size = max_buffer_size &
1103  + pmc_mpi_pack_size_run_part_opt(run_part_opt)
1104  max_buffer_size = max_buffer_size &
1105  + pmc_mpi_pack_size_real(n_part)
1106  max_buffer_size = max_buffer_size &
1107  + pmc_mpi_pack_size_gas_data(gas_data)
1108  max_buffer_size = max_buffer_size &
1109  + pmc_mpi_pack_size_gas_state(gas_state_init)
1110  max_buffer_size = max_buffer_size &
1111  + pmc_mpi_pack_size_aero_data(aero_data)
1112  max_buffer_size = max_buffer_size &
1113  + pmc_mpi_pack_size_aero_dist(aero_dist_init)
1114  max_buffer_size = max_buffer_size &
1115  + pmc_mpi_pack_size_scenario(scenario)
1116  max_buffer_size = max_buffer_size &
1117  + pmc_mpi_pack_size_env_state(env_state_init)
1118  max_buffer_size = max_buffer_size &
1119  + pmc_mpi_pack_size_integer(rand_init)
1120  max_buffer_size = max_buffer_size &
1121  + pmc_mpi_pack_size_logical(do_restart)
1122  max_buffer_size = max_buffer_size &
1123  + pmc_mpi_pack_size_logical(do_init_equilibrate)
1124  max_buffer_size = max_buffer_size &
1125  + pmc_mpi_pack_size_aero_state(aero_state_init)
1126 
1127  allocate(buffer(max_buffer_size))
1128 
1129  position = 0
1130  call pmc_mpi_pack_run_part_opt(buffer, position, run_part_opt)
1131  call pmc_mpi_pack_real(buffer, position, n_part)
1132  call pmc_mpi_pack_gas_data(buffer, position, gas_data)
1133  call pmc_mpi_pack_gas_state(buffer, position, gas_state_init)
1134  call pmc_mpi_pack_aero_data(buffer, position, aero_data)
1135  call pmc_mpi_pack_aero_dist(buffer, position, aero_dist_init)
1136  call pmc_mpi_pack_scenario(buffer, position, scenario)
1137  call pmc_mpi_pack_env_state(buffer, position, env_state_init)
1138  call pmc_mpi_pack_integer(buffer, position, rand_init)
1139  call pmc_mpi_pack_logical(buffer, position, do_restart)
1140  call pmc_mpi_pack_logical(buffer, position, do_init_equilibrate)
1141  call pmc_mpi_pack_aero_state(buffer, position, aero_state_init)
1142  call assert(181905491, position <= max_buffer_size)
1143  buffer_size = position ! might be less than we allocated
1144  end if
1145 
1146  ! tell everyone the size
1147  call pmc_mpi_bcast_integer(buffer_size)
1148 
1149  if (pmc_mpi_rank() /= 0) then
1150  ! non-root processes allocate space
1151  allocate(buffer(buffer_size))
1152  end if
1153 
1154  ! broadcast data to everyone
1155  call pmc_mpi_bcast_packed(buffer)
1156 
1157  if (pmc_mpi_rank() /= 0) then
1158  ! non-root processes unpack data
1159  position = 0
1160  call pmc_mpi_unpack_run_part_opt(buffer, position, run_part_opt)
1161  call pmc_mpi_unpack_real(buffer, position, n_part)
1162  call pmc_mpi_unpack_gas_data(buffer, position, gas_data)
1163  call pmc_mpi_unpack_gas_state(buffer, position, gas_state_init)
1164  call pmc_mpi_unpack_aero_data(buffer, position, aero_data)
1165  call pmc_mpi_unpack_aero_dist(buffer, position, aero_dist_init)
1166  call pmc_mpi_unpack_scenario(buffer, position, scenario)
1167  call pmc_mpi_unpack_env_state(buffer, position, env_state_init)
1168  call pmc_mpi_unpack_integer(buffer, position, rand_init)
1169  call pmc_mpi_unpack_logical(buffer, position, do_restart)
1170  call pmc_mpi_unpack_logical(buffer, position, do_init_equilibrate)
1171  call pmc_mpi_unpack_aero_state(buffer, position, aero_state_init)
1172  call assert(143770146, position == buffer_size)
1173  end if
1174 
1175  ! free the buffer
1176  deallocate(buffer)
1177 #endif
1178 
1179  end subroutine pmc_mpi_broadcast_run_part
1180 
1181 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1182 
1183 end module pmc_run_part
pmc_run_part::parallel_coag_type_local
integer, parameter parallel_coag_type_local
Type code for local parallel coagulation.
Definition: run_part.F90:44
pmc_util::time_to_string_max_len
character(len=pmc_util_convert_string_len) function time_to_string_max_len(time, max_len)
Convert a time to a string format of maximum length.
Definition: util.F90:929
pmc_run_part::pmc_mpi_pack_run_part_opt
subroutine pmc_mpi_pack_run_part_opt(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: run_part.F90:347
pmc_gas_data::gas_data_t
Constant gas data.
Definition: gas_data.F90:35
pmc_mpi::pmc_mpi_size
integer function pmc_mpi_size()
Returns the total number of processes.
Definition: mpi.F90:134
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_aero_state::aero_state_weight_nummass_source
integer, parameter aero_state_weight_nummass_source
Coupled number/mass weighting by source.
Definition: aero_state.F90:51
pmc_gas_state::gas_state_mix
subroutine gas_state_mix(val)
Average val over all processes.
Definition: gas_state.F90:472
pmc_util::system_clock_time
real(kind=dp) function system_clock_time()
Returns the current system clock time in seconds.
Definition: util.F90:1977
pmc_aero_state::aero_state_weight_flat_specified
integer, parameter aero_state_weight_flat_specified
Flat weighting by specified weight classes.
Definition: aero_state.F90:53
pmc_nucleate
Aerosol nucleation functions.
Definition: nucleate.F90:9
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_aero_data::pmc_mpi_unpack_aero_data
subroutine pmc_mpi_unpack_aero_data(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_data.F90:618
pmc_output::input_state
subroutine input_state(filename, index, time, del_t, i_repeat, uuid, aero_data, aero_state, gas_data, gas_state, env_state)
Read the current state.
Definition: output.F90:498
pmc_mpi::pmc_mpi_pack_size_real
integer function pmc_mpi_pack_size_real(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:385
pmc_spec_file::spec_file_close
subroutine spec_file_close(file)
Close a spec file.
Definition: spec_file.F90:135
pmc_mosaic::mosaic_init
subroutine mosaic_init(env_state, aero_data, del_t, do_optical)
Initialize all MOSAIC data-structures.
Definition: mosaic.F90:39
pmc_util::integer_to_string_max_len
character(len=pmc_util_convert_string_len) function integer_to_string_max_len(val, max_len)
Convert an integer to a string format of maximum length.
Definition: util.F90:853
pmc_scenario
The scenario_t structure and associated subroutines.
Definition: scenario.F90:9
pmc_scenario::scenario_t
Scenario data.
Definition: scenario.F90:54
pmc_gas_data
The gas_data_t structure and associated subroutines.
Definition: gas_data.F90:9
pmc_env_state::pmc_mpi_unpack_env_state
subroutine pmc_mpi_unpack_env_state(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: env_state.F90:468
pmc_run_part::run_part
subroutine run_part(scenario, env_state, aero_data, aero_state, gas_data, gas_state, run_part_opt)
Do a particle-resolved Monte Carlo simulation.
Definition: run_part.F90:127
pmc_scenario::pmc_mpi_pack_scenario
subroutine pmc_mpi_pack_scenario(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: scenario.F90:1063
pmc_coagulation_dist::mc_coag_dist
subroutine mc_coag_dist(coag_kernel_type, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag)
Do coagulation for time del_t.
Definition: coagulation_dist.F90:106
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:135
pmc_mpi::pmc_mpi_rank
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
pmc_mpi::pmc_mpi_pack_logical
subroutine pmc_mpi_pack_logical(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:813
pmc_env_state::env_state_mix
subroutine env_state_mix(val)
Average val over all processes.
Definition: env_state.F90:339
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_coag_kernel
Generic coagulation kernel.
Definition: coag_kernel.F90:9
pmc_condense
Water condensation onto aerosol particles.
Definition: condense.F90:29
pmc_mpi::pmc_mpi_pack_string
subroutine pmc_mpi_pack_string(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:786
pmc_run_part::run_part_opt_t
Options controlling the execution of run_part().
Definition: run_part.F90:49
pmc_output::output_state
subroutine output_state(prefix, output_type, aero_data, aero_state, gas_data, gas_state, env_state, index, time, del_t, i_repeat, record_removals, record_optical, uuid)
Write the current state.
Definition: output.F90:112
pmc_aero_state::aero_state_total_particles_all_procs
integer function aero_state_total_particles_all_procs(aero_state, i_group, i_class)
Returns the total number of particles across all processes.
Definition: aero_state.F90:327
pmc_run_part::run_part_timestep
subroutine run_part_timestep(scenario, env_state, aero_data, aero_state, gas_data, gas_state, run_part_opt, i_time, t_start, last_output_time, last_progress_time, i_output, progress_n_samp, progress_n_coag, progress_n_emit, progress_n_dil_in, progress_n_dil_out, progress_n_nuc)
Do a one particle-resolved Monte Carlo simulation timestep.
Definition: run_part.F90:761
pmc_spec_file::spec_file_read_logical
subroutine spec_file_read_logical(file, name, var)
Read a logical from a spec file that must have a given name.
Definition: spec_file.F90:584
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_scenario::pmc_mpi_unpack_scenario
subroutine pmc_mpi_unpack_scenario(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: scenario.F90:1122
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_mpi::pmc_mpi_pack_real
subroutine pmc_mpi_pack_real(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:761
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:104
pmc_scenario::scenario_update_gas_state
subroutine scenario_update_gas_state(scenario, delta_t, env_state, old_env_state, gas_data, gas_state)
Do gas emissions and background dilution.
Definition: scenario.F90:209
pmc_scenario::scenario_update_aero_state
subroutine scenario_update_aero_state(scenario, delta_t, env_state, old_env_state, aero_data, aero_state, n_emit, n_dil_in, n_dil_out, allow_doubling, allow_halving)
Do emissions and background dilution for a particle aerosol distribution.
Definition: scenario.F90:262
pmc_aero_state
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
pmc_mpi::pmc_mpi_unpack_string
subroutine pmc_mpi_unpack_string(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1214
pmc_scenario::scenario_update_env_state
subroutine scenario_update_env_state(scenario, env_state, time)
Update time-dependent contents of the environment. scenario_init_env_state() should have been called ...
Definition: scenario.F90:134
pmc_coag_kernel::coag_kernel_type_additive
integer, parameter coag_kernel_type_additive
Type code for an additive kernel.
Definition: coag_kernel.F90:33
pmc_mpi::pmc_mpi_pack_size_logical
integer function pmc_mpi_pack_size_logical(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:427
pmc_mosaic
Interface to the MOSAIC aerosol and gas phase chemistry code.
Definition: mosaic.F90:9
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_condense::condense_particles
subroutine condense_particles(aero_state, aero_data, env_state_initial, env_state_final, del_t)
Do condensation to all the particles for a given time interval, including updating the environment to...
Definition: condense.F90:149
pmc_mpi::pmc_mpi_bcast_integer
subroutine pmc_mpi_bcast_integer(val)
Broadcast the given value from process 0 to all other processes.
Definition: mpi.F90:288
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:547
pmc_env_state::pmc_mpi_pack_env_state
subroutine pmc_mpi_pack_env_state(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: env_state.F90:421
pmc_run_part::parallel_coag_type_invalid
integer, parameter parallel_coag_type_invalid
Type code for undefined or invalid parallel coagulation method.
Definition: run_part.F90:42
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_util::check_time_multiple
subroutine check_time_multiple(first_name, first_time, second_name, second_time)
Check that the first time interval is close to an integer multiple of the second, and warn if it is n...
Definition: util.F90:377
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:767
pmc_coagulation
Aerosol particle coagulation.
Definition: coagulation.F90:11
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:78
pmc_output::output_type_single
integer, parameter output_type_single
Type code for single output (one file for all processes, written by process 0).
Definition: output.F90:97
pmc_run_part::pmc_mpi_broadcast_run_part
subroutine pmc_mpi_broadcast_run_part(run_part_opt, aero_data, aero_state_init, gas_data, gas_state_init, env_state_init, aero_dist_init, scenario, n_part, rand_init, do_init_equilibrate, do_restart)
Read the specification for a run_part simulation from a spec file.
Definition: run_part.F90:1060
pmc_aero_state::aero_state_total_particles
integer function aero_state_total_particles(aero_state, i_group, i_class)
Returns the total number of particles in an aerosol distribution.
Definition: aero_state.F90:288
pmc_run_part::parallel_coag_type_dist
integer, parameter parallel_coag_type_dist
Type code for distributed parallel coagulation.
Definition: run_part.F90:46
pmc_nucleate::nucleate
subroutine nucleate(nucleate_type, nucleate_source, nucleate_weight_class, env_state, gas_data, aero_data, aero_state, gas_state, del_t, allow_doubling, allow_halving)
Do nucleation of the type given by the first argument.
Definition: nucleate.F90:34
pmc_spec_file::spec_file_open
subroutine spec_file_open(filename, file)
Open a spec file for reading.
Definition: spec_file.F90:112
pmc_scenario::pmc_mpi_pack_size_scenario
integer function pmc_mpi_pack_size_scenario(val)
Determines the number of bytes required to pack the given value.
Definition: scenario.F90:1007
pmc_mosaic::mosaic_cleanup
subroutine mosaic_cleanup()
Clean-up after running MOSAIC, deallocating memory.
Definition: mosaic.F90:131
pmc_photolysis
The photolysis_t type and related functions.
Definition: photolysis.F90:9
pmc_mosaic::mosaic_optical_wavelengths
subroutine mosaic_optical_wavelengths(aero_data)
Fetch the wavelengths that optical cross-sections are calculated at in MOSAIC. If not using multiple ...
Definition: mosaic.F90:679
pmc_run_part
Monte Carlo simulation.
Definition: run_part.F90:9
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_coagulation::mc_coag
subroutine mc_coag(coag_kernel_type, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag)
Do coagulation for time del_t.
Definition: coagulation.F90:45
pmc_mpi::pmc_mpi_bcast_packed
subroutine pmc_mpi_bcast_packed(val)
Broadcast the given value from process 0 to all other processes.
Definition: mpi.F90:326
pmc_tchem_interface
An interface between PartMC and TChem.
Definition: tchem_interface.F90:9
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
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:1139
pmc_mpi::pmc_mpi_unpack_logical
subroutine pmc_mpi_unpack_logical(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1242
pmc_output
Write data in NetCDF format.
Definition: output.F90:68
pmc_constants
Physical constants.
Definition: constants.F90:9
pmc_aero_state::pmc_mpi_unpack_aero_state
subroutine pmc_mpi_unpack_aero_state(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_state.F90:2391
pmc_coag_kernel::coag_kernel_type_invalid
integer, parameter coag_kernel_type_invalid
Type code for an undefined or invalid kernel.
Definition: coag_kernel.F90:29
pmc_run_part::run_part_timeblock
subroutine run_part_timeblock(scenario, env_state, aero_data, aero_state, gas_data, gas_state, run_part_opt, i_cur, i_next, t_start, last_output_time, last_progress_time, i_output, progress_n_samp, progress_n_coag, progress_n_emit, progress_n_dil_in, progress_n_dil_out, progress_n_nuc)
Do a number of time steps of particle-reoslved Monte Carlo simulation.
Definition: run_part.F90:989
pmc_aero_state::pmc_mpi_pack_aero_state
subroutine pmc_mpi_pack_aero_state(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_state.F90:2365
pmc_run_part::spec_file_read_run_part
subroutine spec_file_read_run_part(file, run_part_opt, aero_data, aero_state_init, gas_data, gas_state_init, env_state_init, aero_dist_init, scenario, n_part, rand_init, do_init_equilibrate, do_restart)
Read the specification for a run_part simulation from a spec file.
Definition: run_part.F90:461
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_aero_state::pmc_mpi_pack_size_aero_state
integer function pmc_mpi_pack_size_aero_state(val)
Determines the number of bytes required to pack the given value.
Definition: aero_state.F90:2346
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_util::check_event
subroutine check_event(time, timestep, interval, last_time, do_event)
Computes whether an event is scheduled to take place.
Definition: util.F90:408
pmc_aero_state::aero_state_rebalance
subroutine aero_state_rebalance(aero_state, aero_data, allow_doubling, allow_halving, initial_state_warning)
Double or halve the particle population in each weight group to maintain close to n_part_ideal partic...
Definition: aero_state.F90:1767
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:524
pmc_spec_file::spec_file_die_msg
subroutine spec_file_die_msg(code, file, msg)
Exit with an error message containing filename and line number.
Definition: spec_file.F90:74
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_run_part::pmc_mpi_unpack_run_part_opt
subroutine pmc_mpi_unpack_run_part_opt(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: run_part.F90:401
pmc_camp_interface
An interface between PartMC and the CAMP.
Definition: camp_interface.F90:9
pmc_spec_file::spec_file_read_string
subroutine spec_file_read_string(file, name, var)
Read a string from a spec file that must have a given name.
Definition: spec_file.F90:605
pmc_aero_state::aero_state_mix
subroutine aero_state_mix(aero_state, del_t, mix_timescale, aero_data, specify_prob_transfer)
Mix the aero_states between all processes. Currently uses a simple all-to-all diffusion.
Definition: aero_state.F90:1903
pmc_mosaic::mosaic_aero_optical_init
subroutine mosaic_aero_optical_init(env_state, aero_data, aero_state, gas_data, gas_state)
Compute the optical properties of each aerosol particle for initial timestep.
Definition: mosaic.F90:564
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_mosaic::mosaic_timestep
subroutine mosaic_timestep(env_state, aero_data, aero_state, gas_data, gas_state, do_optical, uuid)
Do one timestep with MOSAIC.
Definition: mosaic.F90:413
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:711
pmc_run_part::print_part_progress
subroutine print_part_progress(i_repeat, t_sim_elapsed, n_part, n_coag, n_emit, n_dil_in, n_dil_out, n_nuc, t_wall_elapsed, t_wall_remain)
Print the current simulation progress to the screen.
Definition: run_part.F90:259
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_mpi::pmc_mpi_reduce_sum_integer
subroutine pmc_mpi_reduce_sum_integer(val, val_sum)
Computes the sum of val across all processes, storing the result in val_sum on the root process.
Definition: mpi.F90:1663
pmc_aero_data::pmc_mpi_pack_size_aero_data
integer function pmc_mpi_pack_size_aero_data(val)
Determines the number of bytes required to pack the given value.
Definition: aero_data.F90:562
pmc_env_state::pmc_mpi_pack_size_env_state
integer function pmc_mpi_pack_size_env_state(val)
Determines the number of bytes required to pack the given value.
Definition: env_state.F90:384
pmc_mpi::pmc_mpi_unpack_real
subroutine pmc_mpi_unpack_real(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1189
pmc_run_part::pmc_mpi_pack_size_run_part_opt
integer function pmc_mpi_pack_size_run_part_opt(val)
Determines the number of bytes required to pack the given value.
Definition: run_part.F90:303
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:511
pmc_nucleate::nucleate_type_invalid
integer, parameter nucleate_type_invalid
Type code for unknown or invalid nucleation type.
Definition: nucleate.F90:18
pmc_aero_state::aero_state_t
The current collection of aerosol particles.
Definition: aero_state.F90:69
pmc_coagulation_dist
Parallel aerosol particle coagulation with MPI.
Definition: coagulation_dist.F90:9
pmc_mosaic::mosaic_support
logical function mosaic_support()
Whether MOSAIC support is compiled in.
Definition: mosaic.F90:26
pmc_mpi::pmc_mpi_pack_size_string
integer function pmc_mpi_pack_size_string(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:405
pmc_aero_data::pmc_mpi_pack_aero_data
subroutine pmc_mpi_pack_aero_data(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_data.F90:585