PartMC  2.6.1
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
28 #ifdef PMC_USE_CAMP
29  use camp_camp_core
30  use camp_camp_state
31  use camp_env_state, only: camp_env_state_t => env_state_t
32 #endif
33 #ifdef PMC_USE_SUNDIALS
34  use pmc_condense
35 #endif
36 #ifdef PMC_USE_MPI
37  use mpi
38 #endif
39 
40  !> Type code for undefined or invalid parallel coagulation method.
41  integer, parameter :: parallel_coag_type_invalid = 0
42  !> Type code for local parallel coagulation.
43  integer, parameter :: parallel_coag_type_local = 1
44  !> Type code for distributed parallel coagulation.
45  integer, parameter :: parallel_coag_type_dist = 2
46 
47  !> Options controlling the execution of run_part().
49  !> Final time (s).
50  real(kind=dp) :: t_max
51  !> Output interval (0 disables) (s).
52  real(kind=dp) :: t_output
53  !> Progress interval (0 disables) (s).
54  real(kind=dp) :: t_progress
55  !> Timestep for coagulation.
56  real(kind=dp) :: del_t
57  !> Prefix for output files.
58  character(len=300) :: output_prefix
59  !> Type of coagulation kernel.
60  integer :: coag_kernel_type
61  !> Type of nucleation.
62  integer :: nucleate_type
63  !> Source of nucleation.
64  integer :: nucleate_source
65  !> Whether to do coagulation.
66  logical :: do_coagulation
67  !> Whether to do nucleation.
68  logical :: do_nucleation
69  !> Allow doubling if needed.
70  logical :: allow_doubling
71  !> Allow halving if needed.
72  logical :: allow_halving
73  !> Whether to do condensation.
74  logical :: do_condensation
75  !> Whether to do MOSAIC.
76  logical :: do_mosaic
77  !> Whether to compute optical properties.
78  logical :: do_optical
79  !> Whether to have explicitly selected weighting.
80  logical :: do_select_weighting
81  !> Type of particle weighting scheme.
82  integer :: weighting_type
83  !> Weighting exponent for power weighting scheme.
84  real(kind=dp) :: weighting_exponent
85  !> Repeat number of run.
86  integer :: i_repeat
87  !> Total number of repeats.
88  integer :: n_repeat
89  !> Cpu_time() of start.
90  real(kind=dp) :: t_wall_start
91  !> Whether to record particle removal information.
92  logical :: record_removals
93  !> Whether to run in parallel.
94  logical :: do_parallel
95  !> Parallel output type.
96  integer :: output_type
97  !> Mixing timescale between processes (s).
98  real(kind=dp) :: mix_timescale
99  !> Whether to average gases each timestep.
100  logical :: gas_average
101  !> Whether to average environment each timestep.
102  logical :: env_average
103  !> Parallel coagulation method type.
104  integer :: parallel_coag_type
105  !> Whether to run CAMP
106  logical :: do_camp_chem
107  !> UUID for this simulation.
108  character(len=PMC_UUID_LEN) :: uuid
109  end type run_part_opt_t
110 
111 contains
112 
113 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
114 
115  !> Do a particle-resolved Monte Carlo simulation.
116 #ifdef PMC_USE_CAMP
117  subroutine run_part(scenario, env_state, aero_data, aero_state, gas_data, &
118  gas_state, run_part_opt, camp_core, photolysis)
119 #else
120  subroutine run_part(scenario, env_state, aero_data, aero_state, gas_data, &
121  gas_state, run_part_opt)
122 #endif
123  !> Environment state.
124  type(scenario_t), intent(in) :: scenario
125  !> Environment state.
126  type(env_state_t), intent(inout) :: env_state
127  !> Aerosol data.
128  type(aero_data_t), intent(in) :: aero_data
129  !> Aerosol state.
130  type(aero_state_t), intent(inout) :: aero_state
131  !> Gas data.
132  type(gas_data_t), intent(in) :: gas_data
133  !> Gas state.
134  type(gas_state_t), intent(inout) :: gas_state
135  !> Monte Carlo options.
136  type(run_part_opt_t), intent(in) :: run_part_opt
137 #ifdef PMC_USE_CAMP
138  !> CAMP chemistry core
139  type(camp_core_t), pointer, intent(inout), optional :: camp_core
140  !> Photolysis calculator
141  type(photolysis_t), pointer, intent(inout), optional :: photolysis
142 
143  type(camp_state_t), pointer :: camp_state
144  type(camp_state_t), pointer :: camp_pre_aero_state, camp_post_aero_state
145 #endif
146 
147  real(kind=dp) :: time, pre_time, pre_del_t, prop_done
148  real(kind=dp) :: last_output_time, last_progress_time
149  integer :: rank, n_proc, pre_index, ncid
150  integer :: pre_i_repeat
151  integer :: n_samp, n_coag, n_emit, n_dil_in, n_dil_out, n_nuc
152  integer :: progress_n_samp, progress_n_coag
153  integer :: progress_n_emit, progress_n_dil_in, progress_n_dil_out
154  integer :: progress_n_nuc, n_part_before
155  integer :: global_n_part, global_n_samp, global_n_coag
156  integer :: global_n_emit, global_n_dil_in, global_n_dil_out
157  integer :: global_n_nuc
158  logical :: do_output, do_state, do_state_netcdf, do_progress, did_coag
159  real(kind=dp) :: t_start, t_wall_now, t_wall_elapsed, t_wall_remain
160  type(env_state_t) :: old_env_state
161  integer :: n_time, i_time, i_time_start, pre_i_time
162  integer :: i_state, i_state_netcdf, i_output
163 #ifdef PMC_USE_CAMP
164  type(camp_env_state_t) :: camp_env_state
165 #endif
166 
167  rank = pmc_mpi_rank()
168  n_proc = pmc_mpi_size()
169 
170  i_time = 0
171  i_output = 1
172  i_state = 1
173  i_state_netcdf = 1
174  time = 0d0
175  progress_n_samp = 0
176  progress_n_coag = 0
177  progress_n_emit = 0
178  progress_n_dil_in = 0
179  progress_n_dil_out = 0
180  progress_n_nuc = 0
181 
182  call check_time_multiple("t_max", run_part_opt%t_max, &
183  "del_t", run_part_opt%del_t)
184  call check_time_multiple("t_output", run_part_opt%t_output, &
185  "del_t", run_part_opt%del_t)
186  call check_time_multiple("t_progress", run_part_opt%t_progress, &
187  "del_t", run_part_opt%del_t)
188 #ifdef PMC_USE_CAMP
189  if (run_part_opt%do_camp_chem) then
190  camp_env_state%temp = env_state%temp
191  camp_env_state%pressure = env_state%pressure
192  camp_state => camp_core%new_state_one_cell(camp_env_state)
193  camp_pre_aero_state => camp_core%new_state_one_cell(camp_env_state)
194  camp_post_aero_state => camp_core%new_state_one_cell(camp_env_state)
195  end if
196 #endif
197 
198  if (run_part_opt%do_mosaic) then
199  call mosaic_init(env_state, aero_data, run_part_opt%del_t, &
200  run_part_opt%do_optical)
201  if (run_part_opt%do_optical) then
202  call mosaic_aero_optical_init(env_state, aero_data, &
203  aero_state, gas_data, gas_state)
204  end if
205  end if
206 
207  if (run_part_opt%t_output > 0d0) then
208  call output_state(run_part_opt%output_prefix, &
209  run_part_opt%output_type, aero_data, aero_state, gas_data, &
210  gas_state, env_state, i_state, time, run_part_opt%del_t, &
211  run_part_opt%i_repeat, run_part_opt%record_removals, &
212  run_part_opt%do_optical, run_part_opt%uuid)
213  call aero_info_array_zero(aero_state%aero_info_array)
214  end if
215 
216  call aero_state_rebalance(aero_state, aero_data, &
217  run_part_opt%allow_doubling, &
218  run_part_opt%allow_halving, initial_state_warning=.true.)
219 
220  t_start = env_state%elapsed_time
221  last_output_time = time
222  last_progress_time = time
223  n_time = nint(run_part_opt%t_max / run_part_opt%del_t)
224  i_time_start = nint(time / run_part_opt%del_t) + 1
225 
226  global_n_part = aero_state_total_particles_all_procs(aero_state)
227  if (rank == 0) then
228  ! progress only printed from root process
229  if (run_part_opt%i_repeat == 1) then
230  t_wall_elapsed = 0d0
231  t_wall_remain = 0d0
232  else
233  call cpu_time(t_wall_now)
234  prop_done = real(run_part_opt%i_repeat - 1, kind=dp) &
235  / real(run_part_opt%n_repeat, kind=dp)
236  t_wall_elapsed = t_wall_now - run_part_opt%t_wall_start
237  t_wall_remain = (1d0 - prop_done) / prop_done &
238  * t_wall_elapsed
239  end if
240  call print_part_progress(run_part_opt%i_repeat, time, &
241  global_n_part, 0, 0, 0, 0, 0, t_wall_elapsed, t_wall_remain)
242  end if
243 
244  do i_time = i_time_start,n_time
245 
246  time = real(i_time, kind=dp) * run_part_opt%del_t
247 
248  old_env_state = env_state
249  call scenario_update_env_state(scenario, env_state, time + t_start)
250 
251  if (run_part_opt%do_nucleation) then
252  n_part_before = aero_state_total_particles(aero_state)
253  call nucleate(run_part_opt%nucleate_type, &
254  run_part_opt%nucleate_source, env_state, gas_data, aero_data, &
255  aero_state, gas_state, run_part_opt%del_t, &
256  run_part_opt%allow_doubling, run_part_opt%allow_halving)
257  n_nuc = aero_state_total_particles(aero_state) &
258  - n_part_before
259  progress_n_nuc = progress_n_nuc + n_nuc
260  end if
261 
262  if (run_part_opt%do_coagulation) then
263  if (run_part_opt%parallel_coag_type &
264  == parallel_coag_type_local) then
265  call mc_coag(run_part_opt%coag_kernel_type, env_state, &
266  aero_data, aero_state, run_part_opt%del_t, n_samp, n_coag)
267  elseif (run_part_opt%parallel_coag_type &
268  == parallel_coag_type_dist) then
269  call mc_coag_dist(run_part_opt%coag_kernel_type, env_state, &
270  aero_data, aero_state, run_part_opt%del_t, n_samp, n_coag)
271  else
272  call die_msg(323011762, "unknown parallel coagulation type: " &
273  // trim(integer_to_string(run_part_opt%parallel_coag_type)))
274  end if
275  progress_n_samp = progress_n_samp + n_samp
276  progress_n_coag = progress_n_coag + n_coag
277  end if
278 
279 #ifdef PMC_USE_SUNDIALS
280  if (run_part_opt%do_condensation) then
281  call condense_particles(aero_state, aero_data, old_env_state, &
282  env_state, run_part_opt%del_t)
283  end if
284 #endif
285 
286  call scenario_update_gas_state(scenario, run_part_opt%del_t, &
287  env_state, old_env_state, gas_data, gas_state)
288  call scenario_update_aero_state(scenario, run_part_opt%del_t, &
289  env_state, old_env_state, aero_data, aero_state, n_emit, &
290  n_dil_in, n_dil_out, run_part_opt%allow_doubling, &
291  run_part_opt%allow_halving)
292  progress_n_emit = progress_n_emit + n_emit
293  progress_n_dil_in = progress_n_dil_in + n_dil_in
294  progress_n_dil_out = progress_n_dil_out + n_dil_out
295 
296  if (run_part_opt%do_camp_chem) then
297 #ifdef PMC_USE_CAMP
298  call pmc_camp_interface_solve(camp_core, camp_state, &
299  camp_pre_aero_state, camp_post_aero_state, env_state, &
300  aero_data, aero_state, gas_data, gas_state, photolysis, &
301  run_part_opt%del_t)
302 #endif
303  end if
304 
305  if (run_part_opt%do_mosaic) then
306  call mosaic_timestep(env_state, aero_data, aero_state, gas_data, &
307  gas_state, run_part_opt%do_optical)
308  end if
309 
310  if (run_part_opt%mix_timescale > 0d0) then
311  call aero_state_mix(aero_state, run_part_opt%del_t, &
312  run_part_opt%mix_timescale, aero_data)
313  end if
314  if (run_part_opt%gas_average) then
315  call gas_state_mix(gas_state)
316  end if
317  if (run_part_opt%gas_average) then
318  call env_state_mix(env_state)
319  end if
320 
321  call aero_state_rebalance(aero_state, aero_data, &
322  run_part_opt%allow_doubling, &
323  run_part_opt%allow_halving, initial_state_warning=.false.)
324 
325  ! DEBUG: enable to check consistency
326  ! call aero_state_check(aero_state, aero_data)
327  ! DEBUG: end
328 
329  if (run_part_opt%t_output > 0d0) then
330  call check_event(time, run_part_opt%del_t, run_part_opt%t_output, &
331  last_output_time, do_output)
332  if (do_output) then
333  i_output = i_output + 1
334  call output_state(run_part_opt%output_prefix, &
335  run_part_opt%output_type, aero_data, aero_state, gas_data, &
336  gas_state, env_state, i_output, time, run_part_opt%del_t, &
337  run_part_opt%i_repeat, run_part_opt%record_removals, &
338  run_part_opt%do_optical, run_part_opt%uuid)
339  call aero_info_array_zero(aero_state%aero_info_array)
340  end if
341  end if
342 
343  if (.not. run_part_opt%record_removals) then
344  ! If we are not recording removals then we can zero them as
345  ! often as possible to minimize the cost of maintaining
346  ! them.
347  call aero_info_array_zero(aero_state%aero_info_array)
348  end if
349 
350  if (run_part_opt%t_progress > 0d0) then
351  call check_event(time, run_part_opt%del_t, &
352  run_part_opt%t_progress, last_progress_time, do_progress)
353  if (do_progress) then
354  global_n_part = aero_state_total_particles_all_procs(aero_state)
355  call pmc_mpi_reduce_sum_integer(progress_n_samp, global_n_samp)
356  call pmc_mpi_reduce_sum_integer(progress_n_coag, global_n_coag)
357  call pmc_mpi_reduce_sum_integer(progress_n_emit, global_n_emit)
358  call pmc_mpi_reduce_sum_integer(progress_n_dil_in, &
359  global_n_dil_in)
360  call pmc_mpi_reduce_sum_integer(progress_n_dil_out, &
361  global_n_dil_out)
362  call pmc_mpi_reduce_sum_integer(progress_n_nuc, global_n_nuc)
363  if (rank == 0) then
364  ! progress only printed from root process
365  call cpu_time(t_wall_now)
366  prop_done = (real(run_part_opt%i_repeat - 1, kind=dp) &
367  + time / run_part_opt%t_max) &
368  / real(run_part_opt%n_repeat, kind=dp)
369  t_wall_elapsed = t_wall_now - run_part_opt%t_wall_start
370  t_wall_remain = (1d0 - prop_done) / prop_done &
371  * t_wall_elapsed
372  call print_part_progress(run_part_opt%i_repeat, time, &
373  global_n_part, global_n_coag, global_n_emit, &
374  global_n_dil_in, global_n_dil_out, global_n_nuc, &
375  t_wall_elapsed, t_wall_remain)
376  end if
377  ! reset counters so they show information since last
378  ! progress display
379  progress_n_samp = 0
380  progress_n_coag = 0
381  progress_n_emit = 0
382  progress_n_dil_in = 0
383  progress_n_dil_out = 0
384  progress_n_nuc = 0
385  end if
386  end if
387 
388  end do
389 
390  if (run_part_opt%do_mosaic) then
391  call mosaic_cleanup()
392  end if
393 
394  end subroutine run_part
395 
396 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
397 
398  !> Print the current simulation progress to the screen.
399  subroutine print_part_progress(i_repeat, t_sim_elapsed, n_part, n_coag, &
400  n_emit, n_dil_in, n_dil_out, n_nuc, t_wall_elapsed, t_wall_remain)
401 
402  !> Repeat number of simulation.
403  integer, intent(in) :: i_repeat
404  !> Elapsed simulation time (s).
405  real(kind=dp), intent(in) :: t_sim_elapsed
406  !> Number of particles.
407  integer, intent(in) :: n_part
408  !> Number of coagulated particles since last progress printing.
409  integer, intent(in) :: n_coag
410  !> Number of emitted particles since last progress printing.
411  integer, intent(in) :: n_emit
412  !> Number of diluted-in particles since last progress printing.
413  integer, intent(in) :: n_dil_in
414  !> Number of diluted-out particles since last progress printing.
415  integer, intent(in) :: n_dil_out
416  !> Number of nucleated particles since last progress printing.
417  integer, intent(in) :: n_nuc
418  !> Elapsed wall time (s).
419  real(kind=dp), intent(in) :: t_wall_elapsed
420  !> Estimated remaining wall time (s).
421  real(kind=dp), intent(in) :: t_wall_remain
422 
423  write(*,'(a6,a1,a6,a1,a7,a1,a7,a1,a7,a1,a8,a1,a9,a1,a7,a1,a6,a1,a6)') &
424  "repeat", " ", "t_sim", " ", "n_part", " ", "n_coag", " ", &
425  "n_emit", " ", "n_dil_in", " ", "n_dil_out", " ", "n_nuc", " ", &
426  "t_wall", " ", "t_rem"
427  write(*,'(a6,a1,a6,a1,a7,a1,a7,a1,a7,a1,a8,a1,a9,a1,a7,a1,a6,a1,a6)') &
428  trim(integer_to_string_max_len(i_repeat, 6)), " ", &
429  trim(time_to_string_max_len(t_sim_elapsed, 6)), " ", &
430  trim(integer_to_string_max_len(n_part, 7)), " ", &
431  trim(integer_to_string_max_len(n_coag, 7)), " ", &
432  trim(integer_to_string_max_len(n_emit, 7)), " ", &
433  trim(integer_to_string_max_len(n_dil_in, 7)), " ", &
434  trim(integer_to_string_max_len(n_dil_out, 7)), " ", &
435  trim(integer_to_string_max_len(n_nuc, 7)), " ", &
436  trim(time_to_string_max_len(t_wall_elapsed, 6)), " ", &
437  trim(time_to_string_max_len(t_wall_remain, 6))
438 
439  end subroutine print_part_progress
440 
441 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
442 
443  !> Determines the number of bytes required to pack the given value.
444  integer function pmc_mpi_pack_size_run_part_opt(val)
445 
446  !> Value to pack.
447  type(run_part_opt_t), intent(in) :: val
448 
450  pmc_mpi_pack_size_real(val%t_max) &
451  + pmc_mpi_pack_size_real(val%t_output) &
452  + pmc_mpi_pack_size_real(val%t_progress) &
453  + pmc_mpi_pack_size_real(val%del_t) &
454  + pmc_mpi_pack_size_string(val%output_prefix) &
455  + pmc_mpi_pack_size_integer(val%coag_kernel_type) &
456  + pmc_mpi_pack_size_integer(val%nucleate_type) &
457  + pmc_mpi_pack_size_integer(val%nucleate_source) &
458  + pmc_mpi_pack_size_logical(val%do_coagulation) &
459  + pmc_mpi_pack_size_logical(val%do_nucleation) &
460  + pmc_mpi_pack_size_logical(val%allow_doubling) &
461  + pmc_mpi_pack_size_logical(val%allow_halving) &
462  + pmc_mpi_pack_size_logical(val%do_condensation) &
463  + pmc_mpi_pack_size_logical(val%do_mosaic) &
464  + pmc_mpi_pack_size_logical(val%do_optical) &
465  + pmc_mpi_pack_size_logical(val%do_select_weighting) &
466  + pmc_mpi_pack_size_integer(val%weighting_type) &
467  + pmc_mpi_pack_size_real(val%weighting_exponent) &
468  + pmc_mpi_pack_size_integer(val%i_repeat) &
469  + pmc_mpi_pack_size_integer(val%n_repeat) &
470  + pmc_mpi_pack_size_real(val%t_wall_start) &
471  + pmc_mpi_pack_size_logical(val%record_removals) &
472  + pmc_mpi_pack_size_logical(val%do_parallel) &
473  + pmc_mpi_pack_size_integer(val%output_type) &
474  + pmc_mpi_pack_size_real(val%mix_timescale) &
475  + pmc_mpi_pack_size_logical(val%gas_average) &
476  + pmc_mpi_pack_size_logical(val%env_average) &
477  + pmc_mpi_pack_size_integer(val%parallel_coag_type) &
478  + pmc_mpi_pack_size_string(val%uuid)
479 
480  end function pmc_mpi_pack_size_run_part_opt
481 
482 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
483 
484  !> Packs the given value into the buffer, advancing position.
485  subroutine pmc_mpi_pack_run_part_opt(buffer, position, val)
486 
487  !> Memory buffer.
488  character, intent(inout) :: buffer(:)
489  !> Current buffer position.
490  integer, intent(inout) :: position
491  !> Value to pack.
492  type(run_part_opt_t), intent(in) :: val
493 
494 #ifdef PMC_USE_MPI
495  integer :: prev_position
496 
497  prev_position = position
498  call pmc_mpi_pack_real(buffer, position, val%t_max)
499  call pmc_mpi_pack_real(buffer, position, val%t_output)
500  call pmc_mpi_pack_real(buffer, position, val%t_progress)
501  call pmc_mpi_pack_real(buffer, position, val%del_t)
502  call pmc_mpi_pack_string(buffer, position, val%output_prefix)
503  call pmc_mpi_pack_integer(buffer, position, val%coag_kernel_type)
504  call pmc_mpi_pack_integer(buffer, position, val%nucleate_type)
505  call pmc_mpi_pack_integer(buffer, position, val%nucleate_source)
506  call pmc_mpi_pack_logical(buffer, position, val%do_coagulation)
507  call pmc_mpi_pack_logical(buffer, position, val%do_nucleation)
508  call pmc_mpi_pack_logical(buffer, position, val%allow_doubling)
509  call pmc_mpi_pack_logical(buffer, position, val%allow_halving)
510  call pmc_mpi_pack_logical(buffer, position, val%do_condensation)
511  call pmc_mpi_pack_logical(buffer, position, val%do_mosaic)
512  call pmc_mpi_pack_logical(buffer, position, val%do_optical)
513  call pmc_mpi_pack_logical(buffer, position, val%do_select_weighting)
514  call pmc_mpi_pack_integer(buffer, position, val%weighting_type)
515  call pmc_mpi_pack_real(buffer, position, val%weighting_exponent)
516  call pmc_mpi_pack_integer(buffer, position, val%i_repeat)
517  call pmc_mpi_pack_integer(buffer, position, val%n_repeat)
518  call pmc_mpi_pack_real(buffer, position, val%t_wall_start)
519  call pmc_mpi_pack_logical(buffer, position, val%record_removals)
520  call pmc_mpi_pack_logical(buffer, position, val%do_parallel)
521  call pmc_mpi_pack_integer(buffer, position, val%output_type)
522  call pmc_mpi_pack_real(buffer, position, val%mix_timescale)
523  call pmc_mpi_pack_logical(buffer, position, val%gas_average)
524  call pmc_mpi_pack_logical(buffer, position, val%env_average)
525  call pmc_mpi_pack_integer(buffer, position, val%parallel_coag_type)
526  call pmc_mpi_pack_string(buffer, position, val%uuid)
527  call assert(946070052, &
528  position - prev_position <= pmc_mpi_pack_size_run_part_opt(val))
529 #endif
530 
531  end subroutine pmc_mpi_pack_run_part_opt
532 
533 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
534 
535  !> Unpacks the given value from the buffer, advancing position.
536  subroutine pmc_mpi_unpack_run_part_opt(buffer, position, val)
537 
538  !> Memory buffer.
539  character, intent(inout) :: buffer(:)
540  !> Current buffer position.
541  integer, intent(inout) :: position
542  !> Value to pack.
543  type(run_part_opt_t), intent(inout) :: val
544 
545 #ifdef PMC_USE_MPI
546  integer :: prev_position
547 
548  prev_position = position
549  call pmc_mpi_unpack_real(buffer, position, val%t_max)
550  call pmc_mpi_unpack_real(buffer, position, val%t_output)
551  call pmc_mpi_unpack_real(buffer, position, val%t_progress)
552  call pmc_mpi_unpack_real(buffer, position, val%del_t)
553  call pmc_mpi_unpack_string(buffer, position, val%output_prefix)
554  call pmc_mpi_unpack_integer(buffer, position, val%coag_kernel_type)
555  call pmc_mpi_unpack_integer(buffer, position, val%nucleate_type)
556  call pmc_mpi_unpack_integer(buffer, position, val%nucleate_source)
557  call pmc_mpi_unpack_logical(buffer, position, val%do_coagulation)
558  call pmc_mpi_unpack_logical(buffer, position, val%do_nucleation)
559  call pmc_mpi_unpack_logical(buffer, position, val%allow_doubling)
560  call pmc_mpi_unpack_logical(buffer, position, val%allow_halving)
561  call pmc_mpi_unpack_logical(buffer, position, val%do_condensation)
562  call pmc_mpi_unpack_logical(buffer, position, val%do_mosaic)
563  call pmc_mpi_unpack_logical(buffer, position, val%do_optical)
564  call pmc_mpi_unpack_logical(buffer, position, val%do_select_weighting)
565  call pmc_mpi_unpack_integer(buffer, position, val%weighting_type)
566  call pmc_mpi_unpack_real(buffer, position, val%weighting_exponent)
567  call pmc_mpi_unpack_integer(buffer, position, val%i_repeat)
568  call pmc_mpi_unpack_integer(buffer, position, val%n_repeat)
569  call pmc_mpi_unpack_real(buffer, position, val%t_wall_start)
570  call pmc_mpi_unpack_logical(buffer, position, val%record_removals)
571  call pmc_mpi_unpack_logical(buffer, position, val%do_parallel)
572  call pmc_mpi_unpack_integer(buffer, position, val%output_type)
573  call pmc_mpi_unpack_real(buffer, position, val%mix_timescale)
574  call pmc_mpi_unpack_logical(buffer, position, val%gas_average)
575  call pmc_mpi_unpack_logical(buffer, position, val%env_average)
576  call pmc_mpi_unpack_integer(buffer, position, val%parallel_coag_type)
577  call pmc_mpi_unpack_string(buffer, position, val%uuid)
578  call assert(480118362, &
579  position - prev_position <= pmc_mpi_pack_size_run_part_opt(val))
580 #endif
581 
582  end subroutine pmc_mpi_unpack_run_part_opt
583 
584 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
585 
586  !> Read the specification for a parallel coagulation type from a spec file.
587  subroutine spec_file_read_parallel_coag_type(file, parallel_coag_type)
588 
589  !> Spec file.
590  type(spec_file_t), intent(inout) :: file
591  !> Kernel type.
592  integer, intent(out) :: parallel_coag_type
593 
594  character(len=SPEC_LINE_MAX_VAR_LEN) :: parallel_coag_type_name
595 
596  !> \page input_format_parallel_coag Input File Format: Parallel Coagulation Type
597  !!
598  !! The output type is specified by the parameter:
599  !! - \b parallel_coag (string): type of parallel coagulation ---
600  !! must be one of: \c local for only within-process
601  !! coagulation or \c dist to have all processes perform
602  !! coagulation globally, requesting particles from other
603  !! processes as needed
604  !!
605  !! See also:
606  !! - \ref spec_file_format --- the input file text format
607 
608  call spec_file_read_string(file, 'parallel_coag', &
609  parallel_coag_type_name)
610  if (trim(parallel_coag_type_name) == 'local') then
611  parallel_coag_type = parallel_coag_type_local
612  elseif (trim(parallel_coag_type_name) == 'dist') then
613  parallel_coag_type = parallel_coag_type_dist
614  else
615  call spec_file_die_msg(494684716, file, &
616  "Unknown parallel coagulation type: " &
617  // trim(parallel_coag_type_name))
618  end if
619 
620  end subroutine spec_file_read_parallel_coag_type
621 
622 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
623 
624 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:43
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:912
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:486
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_gas_state::gas_state_mix
subroutine gas_state_mix(val)
Average val over all processes.
Definition: gas_state.F90:480
pmc_nucleate
Aerosol nucleation functions.
Definition: nucleate.F90:9
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:365
pmc_mosaic::mosaic_init
subroutine mosaic_init(env_state, aero_data, del_t, do_optical)
Initialize all MOSAIC data-structures.
Definition: mosaic.F90:38
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:836
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_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:122
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:134
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:638
pmc_env_state::env_state_mix
subroutine env_state_mix(val)
Average val over all processes.
Definition: env_state.F90:291
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:611
pmc_run_part::run_part_opt_t
Options controlling the execution of run_part().
Definition: run_part.F90:48
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:303
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:586
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
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:259
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:868
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_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:407
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_mosaic::mosaic_timestep
subroutine mosaic_timestep(env_state, aero_data, aero_state, gas_data, gas_state, do_optical)
Do one timestep with MOSAIC.
Definition: mosaic.F90:368
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_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:41
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:376
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:766
pmc_coagulation
Aerosol particle coagulation.
Definition: coagulation.F90:11
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:264
pmc_run_part::parallel_coag_type_dist
integer, parameter parallel_coag_type_dist
Type code for distributed parallel coagulation.
Definition: run_part.F90:45
pmc_mosaic::mosaic_cleanup
subroutine mosaic_cleanup()
Clean-up after running MOSAIC, deallocating memory.
Definition: mosaic.F90:130
pmc_photolysis
The photolysis_t type and related functions.
Definition: photolysis.F90:9
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_nucleate::nucleate
subroutine nucleate(nucleate_type, nucleate_source, 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:33
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:49
pmc_mpi::pmc_mpi_unpack_integer
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:818
pmc_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:896
pmc_output
Write data in NetCDF format.
Definition: output.F90:68
pmc_constants
Physical constants.
Definition: constants.F90:9
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:407
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:1613
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_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:537
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:1743
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:490
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_mpi::pmc_mpi_pack_integer
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:561
pmc_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:401
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:1180
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:843
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:445
pmc_aero_state::aero_state_t
The current collection of aerosol particles.
Definition: aero_state.F90:63
pmc_coagulation_dist
Parallel aerosol particle coagulation with MPI.
Definition: coagulation_dist.F90:9
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:385