source: palm/trunk/SOURCE/time_integration.f90 @ 4576

Last change on this file since 4576 was 4573, checked in by oliver.maas, 4 years ago

pt_surface_heating_rate: calculate pt(0) incrementally by using dt_3d instead of calculating it absolutely by using time_since_reference_point, because time_since_reference_point is set to zero for initializing_actions = 'cyclic_fill', add statement for pt_surface_heating_rate in header.f90

  • Property svn:keywords set to Id
File size: 71.6 KB
Line 
1!> @file time_integration.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: time_integration.f90 4573 2020-06-24 13:08:47Z pavelkrc $
27! calculate pt(0) incrementally by using dt_3d instead of calculating it absolutely
28! by using time_since_reference_point, because time_since_reference_point is set
29! to zero for initializing_actions = 'cyclic_fill'
30!
31! 4565 2020-06-15 08:30:38Z oliver.maas
32! added new surface temperature forcing method for bc_pt_b = 'dirichlet':
33! surface temperature pt(0) can be linearly increased by pt_surface_heating_rate (in K/h)
34!
35! 4564 2020-06-12 14:03:36Z raasch
36! Vertical nesting method of Huq et al. (2019) removed
37!
38! 4521 2020-05-06 11:39:49Z schwenkel
39! Rename variable
40!
41! 4511 2020-04-30 12:20:40Z raasch
42! chemistry decycling replaced by explicit setting of lateral boundary conditions
43!
44! 4508 2020-04-24 13:32:20Z raasch
45! salsa decycling replaced by explicit setting of lateral boundary conditions
46!
47! 4502 2020-04-17 16:14:16Z schwenkel
48! Implementation of ice microphysics
49!
50! 4472 2020-03-24 12:21:00Z Giersch
51! OPENACC COPYIN directive for ddx and ddy added
52!
53! 4466 2020-03-20 16:14:41Z suehring
54! Add advection fluxes to ACC copyin
55!
56! 4457 2020-03-11 14:20:43Z raasch
57! use statement for exchange horiz added
58!
59! 4444 2020-03-05 15:59:50Z raasch
60! bugfix: cpp-directives for serial mode added
61!
62! 4420 2020-02-24 14:13:56Z maronga
63! Added output control for wind turbine model
64!
65! 4403 2020-02-12 13:08:46Z banzhafs
66! Allowing both existing and on-demand emission read modes
67!
68! 4360 2020-01-07 11:25:50Z suehring
69! Bugfix, hour_call_emis uninitialized at first call of time_integration
70!
71! 4346 2019-12-18 11:55:56Z motisi
72! Introduction of wall_flags_total_0, which currently sets bits based on static
73! topography information used in wall_flags_static_0
74!
75! 4329 2019-12-10 15:46:36Z motisi
76! Renamed wall_flags_0 to wall_flags_static_0
77!
78! 4281 2019-10-29 15:15:39Z schwenkel
79! Moved boundary conditions to module interface
80!
81! 4276 2019-10-28 16:03:29Z schwenkel
82! Further modularization of lpm code components
83!
84! 4275 2019-10-28 15:34:55Z schwenkel
85! Move call oft lpm to the end of intermediate timestep loop
86!
87! 4268 2019-10-17 11:29:38Z schwenkel
88! Removing module specific boundary conditions an put them into their modules
89!
90! 4227 2019-09-10 18:04:34Z gronemeier
91! implement new palm_date_time_mod
92!
93! 4226 2019-09-10 17:03:24Z suehring
94! Changes in interface for the offline nesting
95!
96! 4182 2019-08-22 15:20:23Z scharf
97! Corrected "Former revisions" section
98!
99! 4170 2019-08-19 17:12:31Z gronemeier
100! copy diss, diss_p, tdiss_m to GPU
101!
102! 4144 2019-08-06 09:11:47Z raasch
103! relational operators .EQ., .NE., etc. replaced by ==, /=, etc.
104!
105! 4126 2019-07-30 11:09:11Z gronemeier
106! renamed routine to calculate uv exposure
107!
108! 4111 2019-07-22 18:16:57Z suehring
109! advc_flags_1 / advc_flags_2 renamed to advc_flags_m / advc_flags_s
110!
111! 4069 2019-07-01 14:05:51Z Giersch
112! Masked output running index mid has been introduced as a local variable to
113! avoid runtime error (Loop variable has been modified) in time_integration
114!
115! 4064 2019-07-01 05:33:33Z gronemeier
116! Moved call to radiation module out of intermediate time loop
117!
118! 4048 2019-06-21 21:00:21Z knoop
119! Moved production_e_init call into turbulence_closure_mod
120!
121! 4047 2019-06-21 18:58:09Z knoop
122! Added remainings of swap_timelevel upon its dissolution
123!
124! 4043 2019-06-18 16:59:00Z schwenkel
125! Further LPM modularization
126!
127! 4039 2019-06-18 10:32:41Z suehring
128! Rename subroutines in module for diagnostic quantities
129!
130! 4029 2019-06-14 14:04:35Z raasch
131! exchange of ghost points and boundary conditions separated for chemical species and SALSA module,
132! bugfix: decycling of chemistry species after nesting data transfer
133!
134! 4022 2019-06-12 11:52:39Z suehring
135! Call synthetic turbulence generator at last RK3 substep right after boundary
136! conditions are updated in offline nesting in order to assure that
137! perturbations are always imposed
138!
139! 4017 2019-06-06 12:16:46Z schwenkel
140! Mass (volume) flux correction included to ensure global mass conservation for child domains.
141!
142! 3994 2019-05-22 18:08:09Z suehring
143! output of turbulence intensity added
144!
145! 3988 2019-05-22 11:32:37Z kanani
146! Implement steerable output interval for virtual measurements
147!
148! 3968 2019-05-13 11:04:01Z suehring
149! replace nspec_out with n_matched_vars
150!
151! 3929 2019-04-24 12:52:08Z banzhafs
152! Reverse changes back from revision 3878: use chem_boundary_conds instead of
153! chem_boundary_conds_decycle
154!
155!
156! 3885 2019-04-11 11:29:34Z kanani
157! Changes related to global restructuring of location messages and introduction
158! of additional debug messages
159!
160! 3879 2019-04-08 20:25:23Z knoop
161! Moved wtm_forces to module_interface_actions
162!
163! 3872 2019-04-08 15:03:06Z knoop
164! Modifications made for salsa:
165! - Call salsa_emission_update at each time step but do the checks within
166!   salsa_emission_update (i.e. skip_time_do_salsa >= time_since_reference_point
167!   and next_aero_emission_update <= time_since_reference_point ).
168! - Renamed nbins --> nbins_aerosol, ncc_tot --> ncomponents_mass and
169!   ngast --> ngases_salsa and loop indices b, c and sg to ib, ic and ig
170! - Apply nesting for salsa variables
171! - Removed cpu_log calls speciffic for salsa.
172!
173! 3833 2019-03-28 15:04:04Z forkel
174! added USE chem_gasphase_mod, replaced nspec by nspec since fixed compounds are not integrated
175!
176! 3820 2019-03-27 11:53:41Z forkel
177! renamed do_emiss to emissions_anthropogenic (ecc)
178!
179!
180! 3774 2019-03-04 10:52:49Z moh.hefny
181! rephrase if statement to avoid unallocated array in case of
182! nesting_offline is false (crashing during debug mode)
183!
184! 3761 2019-02-25 15:31:42Z raasch $
185! module section re-formatted and openacc required variables moved to separate section,
186! re-formatting to 100 char line width
187!
188! 3745 2019-02-15 18:57:56Z suehring
189! Call indoor model after first timestep
190!
191! 3744 2019-02-15 18:38:58Z suehring
192! - Moved call of bio_calculate_thermal_index_maps from biometeorology module to
193! time_integration to make sure averaged input is updated before calculating.
194!
195! 3739 2019-02-13 08:05:17Z dom_dwd_user
196! Removed everything related to "time_bio_results" as this is never used.
197!
198! 3724 2019-02-06 16:28:23Z kanani
199! Correct double-used log_point_s unit
200!
201! 3719 2019-02-06 13:10:18Z kanani
202! - removed wind_turbine cpu measurement, since same time is measured inside
203!   wtm_forces subroutine as special measures
204! - moved the numerous vnest cpulog to special measures
205! - extended radiation cpulog over entire radiation part,
206!   moved radiation_interactions cpulog to special measures
207! - moved some cpu_log calls to this routine for better overview
208!
209! 3705 2019-01-29 19:56:39Z suehring
210! Data output for virtual measurements added
211!
212! 3704 2019-01-29 19:51:41Z suehring
213! Rename subroutines for surface-data output
214!
215! 3647 2019-01-02 14:10:44Z kanani
216! Bugfix: add time_since_reference_point to IF clause for data_output calls
217! (otherwise skip_time_* values don't come into affect with dt_do* = 0.0).
218! Clean up indoor_model and biometeorology model call.
219!
220! Revision 1.1  1997/08/11 06:19:04  raasch
221! Initial revision
222!
223!
224! Description:
225! ------------
226!> Integration in time of the model equations, statistical analysis and graphic
227!> output
228!------------------------------------------------------------------------------!
229 SUBROUTINE time_integration
230
231
232    USE advec_ws,                                                                                  &
233        ONLY:  ws_statistics
234
235    USE arrays_3d,                                                                                 &
236        ONLY:  diss, diss_p, dzu, e_p, nc_p, ni_p, nr_p, prho, pt, pt_p, pt_init, q, qc_p, qr_p,   &
237               q_init, q_p, qi_p, ref_state, rho_ocean, sa_p, s_p, tend, u, u_p, v, vpt, v_p, w_p
238
239#if defined( __parallel )  &&  ! defined( _OPENACC )
240    USE arrays_3d,                                                                                 &
241        ONLY:  e, nc, ni, nr, qc, qi, qr, s, w
242#endif
243
244    USE biometeorology_mod,                                                                        &
245        ONLY:  bio_calculate_thermal_index_maps, thermal_comfort, bio_calculate_uv_exposure,       &
246               uv_exposure
247
248    USE bulk_cloud_model_mod,                                                                      &
249        ONLY: bulk_cloud_model, calc_liquid_water_content, collision_turbulence,                   &
250              microphysics_ice_phase, microphysics_morrison, microphysics_seifert
251
252    USE calc_mean_profile_mod,                                                                     &
253        ONLY:  calc_mean_profile
254
255    USE chem_emissions_mod,                                                                        &
256        ONLY:  chem_emissions_setup, chem_emissions_update_on_demand
257
258    USE chem_gasphase_mod,                                                                         &
259        ONLY:  nvar
260
261    USE chem_modules,                                                                              &
262        ONLY:  bc_cs_t_val, chem_species, communicator_chem, emissions_anthropogenic,              &
263               emiss_read_legacy_mode, n_matched_vars
264
265    USE control_parameters,                                                                        &
266        ONLY:  advected_distance_x, advected_distance_y, air_chemistry, average_count_3d,          &
267               averaging_interval, averaging_interval_pr, bc_lr_cyc, bc_ns_cyc, bc_pt_t_val,       &
268               bc_q_t_val, biometeorology, call_psolver_at_all_substeps,  child_domain,            &
269               constant_flux_layer, constant_heatflux, create_disturbances,        &
270               dopr_n, constant_diffusion, coupling_mode, coupling_start_time,                     &
271               current_timestep_number, disturbance_created, disturbance_energy_limit, dist_range, &
272               do_sum, dt_3d, dt_averaging_input, dt_averaging_input_pr, dt_coupling,              &
273               dt_data_output_av, dt_disturb, dt_do2d_xy, dt_do2d_xz, dt_do2d_yz, dt_do3d,         &
274               dt_domask,dt_dopts, dt_dopr, dt_dopr_listing, dt_dots, dt_run_control,              &
275               end_time, first_call_lpm, first_call_mas, galilei_transformation, humidity,         &
276               indoor_model, intermediate_timestep_count, intermediate_timestep_count_max,         &
277               land_surface, large_scale_forcing, loop_optimization, lsf_surf, lsf_vert, masks,    &
278               multi_agent_system_end, multi_agent_system_start, nesting_offline, neutral,         &
279               nr_timesteps_this_run, nudging, ocean_mode, passive_scalar, pt_reference,           &
280               pt_slope_offset, pt_surface, pt_surface_heating_rate, pt_surface_initial_change,    &
281               random_heatflux, rans_tke_e, run_coupled, salsa,                                    &
282               simulated_time, simulated_time_chr, skip_time_do2d_xy, skip_time_do2d_xz,           &
283               skip_time_do2d_yz, skip_time_do3d, skip_time_domask, skip_time_dopr,                &
284               skip_time_data_output_av, sloping_surface, stop_dt, surface_output,                 &
285               terminate_coupled, terminate_run, timestep_scheme, time_coupling, time_do2d_xy,     &
286               time_do2d_xz, time_do2d_yz, time_do3d, time_domask, time_dopr, time_dopr_av,        &
287               time_dopr_listing, time_dopts, time_dosp, time_dosp_av, time_dots, time_do_av,      &
288               time_do_sla, time_disturb, time_run_control, time_since_reference_point,            &
289               timestep_count, turbulent_inflow, turbulent_outflow, urban_surface,                 &
290               use_initial_profile_as_reference, use_single_reference_value, u_gtrans, v_gtrans,   &
291               virtual_flight, virtual_measurement, ws_scheme_mom, ws_scheme_sca
292
293#if defined( __parallel )
294    USE control_parameters,                                                                        &
295        ONLY:  rans_mode
296#endif
297
298    USE cpulog,                                                                                    &
299        ONLY:  cpu_log, log_point, log_point_s
300
301    USE diagnostic_output_quantities_mod,                                                          &
302        ONLY:  doq_calculate,                                                                      &
303               timestep_number_at_prev_calc
304
305    USE exchange_horiz_mod,                                                                        &
306        ONLY:  exchange_horiz
307
308    USE flight_mod,                                                                                &
309        ONLY:  flight_measurement
310
311    USE grid_variables,                                                                            &
312        ONLY:  ddx, ddy
313
314    USE indices,                                                                                   &
315        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, nzb, nzt
316
317    USE indoor_model_mod,                                                                          &
318        ONLY:  dt_indoor, im_main_heatcool, time_indoor
319
320    USE interfaces
321
322    USE kinds
323
324    USE land_surface_model_mod,                                                                    &
325        ONLY:  lsm_boundary_condition, lsm_energy_balance, lsm_soil_model, skip_time_do_lsm
326
327    USE lagrangian_particle_model_mod,                                                             &
328        ONLY:  lpm_data_output_ptseries
329
330    USE lsf_nudging_mod,                                                                           &
331        ONLY:  calc_tnudge, ls_forcing_surf, ls_forcing_vert, nudge_ref
332
333    USE module_interface,                                                                          &
334        ONLY:  module_interface_actions, module_interface_swap_timelevel,                          &
335               module_interface_boundary_conditions
336
337    USE multi_agent_system_mod,                                                                    &
338        ONLY:  agents_active, multi_agent_system
339
340    USE nesting_offl_mod,                                                                          &
341        ONLY:  nesting_offl_bc,                                                                    &
342               nesting_offl_geostrophic_wind,                                                      &
343               nesting_offl_input,                                                                 &
344               nesting_offl_interpolation_factor,                                                  &
345               nesting_offl_mass_conservation
346
347    USE netcdf_data_input_mod,                                                                     &
348        ONLY:  chem_emis, chem_emis_att
349
350    USE ocean_mod,                                                                                 &
351        ONLY:  prho_reference
352
353    USE palm_date_time_mod,                                                                        &
354        ONLY:  get_date_time
355
356    USE particle_attributes,                                                                       &
357        ONLY:  particle_advection, particle_advection_start, use_sgs_for_particles, wang_kernel
358
359    USE pegrid
360
361#if defined( __parallel )
362    USE pmc_interface,                                                                             &
363        ONLY:  nested_run, nesting_mode, pmci_boundary_conds, pmci_datatrans, pmci_synchronize,    &
364        pmci_ensure_nest_mass_conservation, pmci_ensure_nest_mass_conservation_vertical,           &
365        pmci_set_swaplevel
366#endif
367
368    USE progress_bar,                                                                              &
369        ONLY:  finish_progress_bar, output_progress_bar
370
371    USE prognostic_equations_mod,                                                                  &
372        ONLY:  prognostic_equations_cache, prognostic_equations_vector
373
374    USE radiation_model_mod,                                                                       &
375        ONLY: dt_radiation, force_radiation_call, radiation, radiation_control,                    &
376              radiation_interaction, radiation_interactions, skip_time_do_radiation, time_radiation
377
378    USE salsa_mod,                                                                                 &
379        ONLY: aerosol_number, aerosol_mass, bc_am_t_val, bc_an_t_val, bc_gt_t_val,                 &
380              communicator_salsa, nbins_aerosol, ncomponents_mass, ngases_salsa,                   &
381              salsa_boundary_conditions, salsa_emission_update, salsa_gas, salsa_gases_from_chem,  &
382              skip_time_do_salsa
383
384    USE spectra_mod,                                                                               &
385        ONLY: average_count_sp, averaging_interval_sp, calc_spectra, dt_dosp, skip_time_dosp
386
387    USE statistics,                                                                                &
388        ONLY:  flow_statistics_called, hom, pr_palm, sums_ls_l
389
390
391    USE surface_layer_fluxes_mod,                                                                  &
392        ONLY:  surface_layer_fluxes
393
394    USE surface_data_output_mod,                                                                   &
395        ONLY:  average_count_surf, averaging_interval_surf, dt_dosurf, dt_dosurf_av,               &
396               surface_data_output, surface_data_output_averaging, skip_time_dosurf,               &
397               skip_time_dosurf_av, time_dosurf, time_dosurf_av
398
399    USE surface_mod,                                                                               &
400        ONLY:  surf_def_h, surf_lsm_h, surf_usm_h
401
402    USE synthetic_turbulence_generator_mod,                                                        &
403        ONLY:  dt_stg_call, dt_stg_adjust, parametrize_inflow_turbulence, stg_adjust, stg_main,    &
404               time_stg_adjust, time_stg_call, use_syn_turb_gen
405
406    USE turbulence_closure_mod,                                                                    &
407        ONLY:  tcm_diffusivities
408
409    USE urban_surface_mod,                                                                         &
410        ONLY:  usm_boundary_condition, usm_material_heat_model, usm_material_model,                &
411               usm_surface_energy_balance, usm_green_heat_model
412
413    USE virtual_measurement_mod,                                                                   &
414        ONLY:  dt_virtual_measurement,                                                             &
415               time_virtual_measurement,                                                           &
416               vm_data_output,                                                                     &
417               vm_sampling,                                                                        &
418               vm_time_start
419
420
421    USE wind_turbine_model_mod,                                                                    &
422        ONLY:  dt_data_output_wtm, time_wtm, wind_turbine, wtm_data_output
423
424#if defined( _OPENACC )
425    USE arrays_3d,                                                                                 &
426        ONLY:  d, dd2zu, ddzu, ddzw,                                                               &
427               diss_l_u,                                                                           &
428               diss_l_v,                                                                           &
429               diss_l_w,                                                                           &
430               diss_s_u,                                                                           &
431               diss_s_v,                                                                           &
432               diss_s_w,                                                                           &
433               drho_air, drho_air_zw, dzw, e,                                                      &
434               flux_l_u,                                                                           &
435               flux_l_v,                                                                           &
436               flux_l_w,                                                                           &
437               flux_s_u,                                                                           &
438               flux_s_v,                                                                           &
439               flux_s_w,                                                                           &
440               heatflux_output_conversion,                                                         &
441               kh, km, momentumflux_output_conversion, nc, ni, nr, p, ptdf_x, ptdf_y, qc, qi, qr, rdf,     &
442               rdf_sc, rho_air, rho_air_zw, s, tdiss_m, te_m, tpt_m, tu_m, tv_m, tw_m, ug, u_init, &
443               u_stokes_zu, vg, v_init, v_stokes_zu, w, zu
444
445    USE control_parameters,                                                                        &
446        ONLY:  tsc
447
448    USE indices,                                                                                   &
449        ONLY:  advc_flags_m, advc_flags_s, nyn, nyng, nys, nysg, nz, nzb_max, wall_flags_total_0
450
451    USE statistics,                                                                                &
452        ONLY:  rmask, statistic_regions, sums_l, sums_l_l, sums_us2_ws_l,                          &
453               sums_wsus_ws_l, sums_vs2_ws_l, sums_wsvs_ws_l, sums_ws2_ws_l, sums_wspts_ws_l,      &
454               sums_wsqs_ws_l, sums_wssas_ws_l, sums_wsqcs_ws_l, sums_wsqrs_ws_l, sums_wsncs_ws_l, &
455               sums_wsnrs_ws_l, sums_wsss_ws_l, weight_substep, sums_salsa_ws_l, sums_wsqis_ws_l,  &
456               sums_wsnis_ws_l
457
458    USE surface_mod,                                                                               &
459        ONLY:  bc_h, enter_surface_arrays, exit_surface_arrays
460#endif
461
462
463    IMPLICIT NONE
464
465    CHARACTER (LEN=9) ::  time_to_string   !<
466
467    INTEGER(iwp) ::  hour                !< hour of current time
468    INTEGER(iwp) ::  hour_call_emis = -1 !< last hour where emission was called
469    INTEGER(iwp) ::  ib                  !< index for aerosol size bins
470    INTEGER(iwp) ::  ic                  !< index for aerosol mass bins
471    INTEGER(iwp) ::  icc                 !< additional index for aerosol mass bins
472    INTEGER(iwp) ::  ig                  !< index for salsa gases
473    INTEGER(iwp) ::  mid                 !< masked output running index
474    INTEGER(iwp) ::  n                   !< loop counter for chemistry species
475
476    REAL(wp) ::  dt_3d_old  !< temporary storage of timestep to be used for
477                            !< steering of run control output interval
478    REAL(wp) ::  time_since_reference_point_save  !< original value of
479                                                  !< time_since_reference_point
480
481
482!
483!-- Copy data from arrays_3d
484!$ACC DATA &
485!$ACC COPY(d(nzb+1:nzt,nys:nyn,nxl:nxr)) &
486!$ACC COPY(diss(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
487!$ACC COPY(e(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
488!$ACC COPY(u(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
489!$ACC COPY(v(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
490!$ACC COPY(w(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
491!$ACC COPY(kh(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
492!$ACC COPY(km(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
493!$ACC COPY(p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
494!$ACC COPY(pt(nzb:nzt+1,nysg:nyng,nxlg:nxrg))
495
496!$ACC DATA &
497!$ACC COPYIN(diss_l_u(0:nz+1,nys:nyn,0), flux_l_u(0:nz+1,nys:nyn,0)) &
498!$ACC COPYIN(diss_l_v(0:nz+1,nys:nyn,0), flux_l_v(0:nz+1,nys:nyn,0)) &
499!$ACC COPYIN(diss_l_w(0:nz+1,nys:nyn,0), flux_l_w(0:nz+1,nys:nyn,0)) &
500!$ACC COPYIN(diss_s_u(0:nz+1,0), flux_s_u(0:nz+1,0)) &
501!$ACC COPYIN(diss_s_v(0:nz+1,0), flux_s_v(0:nz+1,0)) &
502!$ACC COPYIN(diss_s_w(0:nz+1,0), flux_s_w(0:nz+1,0)) &
503!$ACC COPY(diss_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
504!$ACC COPY(e_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
505!$ACC COPY(u_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
506!$ACC COPY(v_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
507!$ACC COPY(w_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
508!$ACC COPY(pt_p(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
509!$ACC COPY(tend(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
510!$ACC COPY(tdiss_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
511!$ACC COPY(te_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
512!$ACC COPY(tu_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
513!$ACC COPY(tv_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
514!$ACC COPY(tw_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
515!$ACC COPY(tpt_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg))
516
517!$ACC DATA &
518!$ACC COPYIN(rho_air(nzb:nzt+1), drho_air(nzb:nzt+1)) &
519!$ACC COPYIN(rho_air_zw(nzb:nzt+1), drho_air_zw(nzb:nzt+1)) &
520!$ACC COPYIN(zu(nzb:nzt+1)) &
521!$ACC COPYIN(dzu(1:nzt+1), dzw(1:nzt+1)) &
522!$ACC COPYIN(ddzu(1:nzt+1), dd2zu(1:nzt)) &
523!$ACC COPYIN(ddzw(1:nzt+1)) &
524!$ACC COPYIN(heatflux_output_conversion(nzb:nzt+1)) &
525!$ACC COPYIN(momentumflux_output_conversion(nzb:nzt+1)) &
526!$ACC COPYIN(rdf(nzb+1:nzt), rdf_sc(nzb+1:nzt)) &
527!$ACC COPYIN(ptdf_x(nxlg:nxrg), ptdf_y(nysg:nyng)) &
528!$ACC COPYIN(ref_state(0:nz+1)) &
529!$ACC COPYIN(u_init(0:nz+1), v_init(0:nz+1)) &
530!$ACC COPYIN(u_stokes_zu(nzb:nzt+1), v_stokes_zu(nzb:nzt+1)) &
531!$ACC COPYIN(pt_init(0:nz+1)) &
532!$ACC COPYIN(ug(0:nz+1), vg(0:nz+1))
533
534!
535!-- Copy data from control_parameters
536!$ACC DATA &
537!$ACC COPYIN(tsc(1:5))
538
539!
540!-- Copy data from grid_variables
541!$ACC DATA &
542!$ACC COPYIN(ddx, ddy)
543
544!
545!-- Copy data from indices
546!$ACC DATA &
547!$ACC COPYIN(advc_flags_m(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
548!$ACC COPYIN(advc_flags_s(nzb:nzt+1,nysg:nyng,nxlg:nxrg)) &
549!$ACC COPYIN(wall_flags_total_0(nzb:nzt+1,nysg:nyng,nxlg:nxrg))
550
551!
552!-- Copy data from surface_mod
553!$ACC DATA &
554!$ACC COPYIN(bc_h(0:1)) &
555!$ACC COPYIN(bc_h(0)%i(1:bc_h(0)%ns)) &
556!$ACC COPYIN(bc_h(0)%j(1:bc_h(0)%ns)) &
557!$ACC COPYIN(bc_h(0)%k(1:bc_h(0)%ns)) &
558!$ACC COPYIN(bc_h(1)%i(1:bc_h(1)%ns)) &
559!$ACC COPYIN(bc_h(1)%j(1:bc_h(1)%ns)) &
560!$ACC COPYIN(bc_h(1)%k(1:bc_h(1)%ns))
561
562!
563!-- Copy data from statistics
564!$ACC DATA &
565!$ACC COPYIN(hom(0:nz+1,1:2,1:4,0)) &
566!$ACC COPYIN(rmask(nysg:nyng,nxlg:nxrg,0:statistic_regions)) &
567!$ACC COPYIN(weight_substep(1:intermediate_timestep_count_max)) &
568!$ACC COPY(sums_l(nzb:nzt+1,1:pr_palm,0)) &
569!$ACC COPY(sums_l_l(nzb:nzt+1,0:statistic_regions,0)) &
570!$ACC COPY(sums_us2_ws_l(nzb:nzt+1,0)) &
571!$ACC COPY(sums_wsus_ws_l(nzb:nzt+1,0)) &
572!$ACC COPY(sums_vs2_ws_l(nzb:nzt+1,0)) &
573!$ACC COPY(sums_wsvs_ws_l(nzb:nzt+1,0)) &
574!$ACC COPY(sums_ws2_ws_l(nzb:nzt+1,0)) &
575!$ACC COPY(sums_wspts_ws_l(nzb:nzt+1,0)) &
576!$ACC COPY(sums_wssas_ws_l(nzb:nzt+1,0)) &
577!$ACC COPY(sums_wsqs_ws_l(nzb:nzt+1,0)) &
578!$ACC COPY(sums_wsqcs_ws_l(nzb:nzt+1,0)) &
579!$ACC COPY(sums_wsqis_ws_l(nzb:nzt+1,0)) &
580!$ACC COPY(sums_wsqrs_ws_l(nzb:nzt+1,0)) &
581!$ACC COPY(sums_wsncs_ws_l(nzb:nzt+1,0)) &
582!$ACC COPY(sums_wsnis_ws_l(nzb:nzt+1,0)) &
583!$ACC COPY(sums_wsnrs_ws_l(nzb:nzt+1,0)) &
584!$ACC COPY(sums_wsss_ws_l(nzb:nzt+1,0)) &
585!$ACC COPY(sums_salsa_ws_l(nzb:nzt+1,0))
586
587!
588!-- Next statement is to avoid compiler warnings about unused variables. Please
589!-- remove in case that you are using them. ddx and ddy need to be defined in
590!-- time_integration because of ACC COPYIN directive.
591    ddx = ddx
592    ddy = ddy
593
594#if defined( _OPENACC )
595    CALL enter_surface_arrays
596#endif
597
598!
599!-- At beginning determine the first time step
600    CALL timestep
601
602#if defined( __parallel )
603!
604!-- Synchronize the timestep in case of nested run.
605    IF ( nested_run )  THEN
606!
607!--    Synchronization by unifying the time step.
608!--    Global minimum of all time-steps is used for all.
609       CALL pmci_synchronize
610    ENDIF
611#endif
612
613!
614!-- Determine and print out the run control quantities before the first time
615!-- step of this run. For the initial run, some statistics (e.g. divergence)
616!-- need to be determined first --> CALL flow_statistics at the beginning of
617!-- run_control
618    CALL run_control
619!
620!-- Data exchange between coupled models in case that a call has been omitted
621!-- at the end of the previous run of a job chain.
622    IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
623!
624!--    In case of model termination initiated by the local model the coupler
625!--    must not be called because this would again cause an MPI hang.
626       DO WHILE ( time_coupling >= dt_coupling  .AND.  terminate_coupled == 0 )
627          CALL surface_coupler
628          time_coupling = time_coupling - dt_coupling
629       ENDDO
630       IF (time_coupling == 0.0_wp  .AND.  time_since_reference_point < dt_coupling )  THEN
631          time_coupling = time_since_reference_point
632       ENDIF
633    ENDIF
634
635    CALL location_message( 'atmosphere (and/or ocean) time-stepping', 'start' )
636
637!
638!-- Start of the time loop
639    DO  WHILE ( simulated_time < end_time  .AND.  .NOT. stop_dt  .AND. .NOT. terminate_run )
640
641       CALL cpu_log( log_point_s(10), 'timesteps', 'start' )
642
643!
644!--    Determine ug, vg and w_subs in dependence on data from external file
645!--    LSF_DATA
646       IF ( large_scale_forcing .AND. lsf_vert )  THEN
647           CALL ls_forcing_vert ( simulated_time )
648           sums_ls_l = 0.0_wp
649       ENDIF
650
651!
652!--    Set pt_init and q_init to the current profiles taken from
653!--    NUDGING_DATA
654       IF ( nudging )  THEN
655           CALL nudge_ref ( simulated_time )
656!
657!--        Store temperature gradient at the top boundary for possible Neumann
658!--        boundary condition
659           bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
660           bc_q_t_val  = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1)
661           IF ( air_chemistry )  THEN
662              DO  n = 1, nvar
663                 bc_cs_t_val = (  chem_species(n)%conc_pr_init(nzt+1)                              &
664                                - chem_species(n)%conc_pr_init(nzt) )                              &
665                               / dzu(nzt+1)
666              ENDDO
667           ENDIF
668           IF ( salsa  .AND.  time_since_reference_point >= skip_time_do_salsa )  THEN
669              DO  ib = 1, nbins_aerosol
670                 bc_an_t_val = ( aerosol_number(ib)%init(nzt+1) - aerosol_number(ib)%init(nzt) ) / &
671                               dzu(nzt+1)
672                 DO  ic = 1, ncomponents_mass
673                    icc = ( ic - 1 ) * nbins_aerosol + ib
674                    bc_am_t_val = ( aerosol_mass(icc)%init(nzt+1) - aerosol_mass(icc)%init(nzt) ) /&
675                                  dzu(nzt+1)
676                 ENDDO
677              ENDDO
678              IF ( .NOT. salsa_gases_from_chem )  THEN
679                 DO  ig = 1, ngases_salsa
680                    bc_gt_t_val = ( salsa_gas(ig)%init(nzt+1) - salsa_gas(ig)%init(nzt) ) /        &
681                                  dzu(nzt+1)
682                 ENDDO
683              ENDIF
684           ENDIF
685       ENDIF
686!
687!--    Input of boundary data.
688       IF ( nesting_offline )  CALL nesting_offl_input
689!
690!--    Execute all other module actions routines
691       CALL module_interface_actions( 'before_timestep' )
692
693!
694!--    Start of intermediate step loop
695       intermediate_timestep_count = 0
696       DO  WHILE ( intermediate_timestep_count < intermediate_timestep_count_max )
697
698          intermediate_timestep_count = intermediate_timestep_count + 1
699
700!
701!--       Set the steering factors for the prognostic equations which depend
702!--       on the timestep scheme
703          CALL timestep_scheme_steering
704
705!
706!--       Calculate those variables needed in the tendency terms which need
707!--       global communication
708          IF ( .NOT. use_single_reference_value  .AND.  .NOT. use_initial_profile_as_reference )   &
709          THEN
710!
711!--          Horizontally averaged profiles to be used as reference state in
712!--          buoyancy terms (WARNING: only the respective last call of
713!--          calc_mean_profile defines the reference state!)
714             IF ( .NOT. neutral )  THEN
715                CALL calc_mean_profile( pt, 4 )
716                ref_state(:)  = hom(:,1,4,0) ! this is used in the buoyancy term
717             ENDIF
718             IF ( ocean_mode )  THEN
719                CALL calc_mean_profile( rho_ocean, 64 )
720                ref_state(:)  = hom(:,1,64,0)
721             ENDIF
722             IF ( humidity )  THEN
723                CALL calc_mean_profile( vpt, 44 )
724                ref_state(:)  = hom(:,1,44,0)
725             ENDIF
726!
727!--          Assure that ref_state does not become zero at any level
728!--          ( might be the case if a vertical level is completely occupied
729!--            with topography ).
730             ref_state = MERGE( MAXVAL(ref_state), ref_state, ref_state == 0.0_wp )
731          ENDIF
732
733          IF ( ( ws_scheme_mom .OR. ws_scheme_sca )  .AND.  intermediate_timestep_count == 1 )     &
734          THEN
735             CALL ws_statistics
736          ENDIF
737!
738!--       In case of nudging calculate current nudging time scale and horizontal
739!--       means of u, v, pt and q
740          IF ( nudging )  THEN
741             CALL calc_tnudge( simulated_time )
742             CALL calc_mean_profile( u, 1 )
743             CALL calc_mean_profile( v, 2 )
744             CALL calc_mean_profile( pt, 4 )
745             CALL calc_mean_profile( q, 41 )
746          ENDIF
747!
748!--       Execute all other module actions routunes
749          CALL module_interface_actions( 'before_prognostic_equations' )
750!
751!--       Solve the prognostic equations. A fast cache optimized version with
752!--       only one single loop is used in case of Piascek-Williams advection
753!--       scheme. NEC vector machines use a different version, because
754!--       in the other versions a good vectorization is prohibited due to
755!--       inlining problems.
756          IF ( loop_optimization == 'cache' )  THEN
757             CALL prognostic_equations_cache
758          ELSEIF ( loop_optimization == 'vector' )  THEN
759             CALL prognostic_equations_vector
760          ENDIF
761
762!
763!--       Movement of agents in multi agent system
764          IF ( agents_active  .AND.  time_since_reference_point >= multi_agent_system_start  .AND. &
765               time_since_reference_point <= multi_agent_system_end  .AND.                         &
766               intermediate_timestep_count == 1 )                                                  &
767          THEN
768             CALL multi_agent_system
769             first_call_mas = .FALSE.
770          ENDIF
771
772!
773!--       Exchange of ghost points (lateral boundary conditions)
774          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'start' )
775
776          CALL exchange_horiz( u_p, nbgp )
777          CALL exchange_horiz( v_p, nbgp )
778          CALL exchange_horiz( w_p, nbgp )
779          CALL exchange_horiz( pt_p, nbgp )
780          IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e_p, nbgp )
781          IF ( rans_tke_e  .OR.  wang_kernel  .OR.  collision_turbulence                           &
782               .OR.  use_sgs_for_particles )  THEN
783             IF ( rans_tke_e )  THEN
784                CALL exchange_horiz( diss_p, nbgp )
785             ELSE
786                CALL exchange_horiz( diss, nbgp )
787             ENDIF
788          ENDIF
789          IF ( ocean_mode )  THEN
790             CALL exchange_horiz( sa_p, nbgp )
791             CALL exchange_horiz( rho_ocean, nbgp )
792             CALL exchange_horiz( prho, nbgp )
793          ENDIF
794          IF ( humidity )  THEN
795             CALL exchange_horiz( q_p, nbgp )
796             IF ( bulk_cloud_model .AND. microphysics_morrison )  THEN
797                CALL exchange_horiz( qc_p, nbgp )
798                CALL exchange_horiz( nc_p, nbgp )
799             ENDIF
800             IF ( bulk_cloud_model .AND. microphysics_seifert )  THEN
801                CALL exchange_horiz( qr_p, nbgp )
802                CALL exchange_horiz( nr_p, nbgp )
803             ENDIF
804             IF ( bulk_cloud_model .AND. microphysics_ice_phase )  THEN
805                CALL exchange_horiz( qi_p, nbgp )
806                CALL exchange_horiz( ni_p, nbgp )
807             ENDIF
808          ENDIF
809          IF ( passive_scalar )  CALL exchange_horiz( s_p, nbgp )
810          IF ( air_chemistry )  THEN
811             DO  n = 1, nvar
812                CALL exchange_horiz( chem_species(n)%conc_p, nbgp,                                 &
813                                     alternative_communicator = communicator_chem )
814             ENDDO
815          ENDIF
816
817          IF ( salsa  .AND.  time_since_reference_point >= skip_time_do_salsa )  THEN
818             DO  ib = 1, nbins_aerosol
819                CALL exchange_horiz( aerosol_number(ib)%conc_p, nbgp,                              &
820                                     alternative_communicator = communicator_salsa )
821                DO  ic = 1, ncomponents_mass
822                   icc = ( ic - 1 ) * nbins_aerosol + ib
823                   CALL exchange_horiz( aerosol_mass(icc)%conc_p, nbgp,                            &
824                                        alternative_communicator = communicator_salsa )
825                ENDDO
826             ENDDO
827             IF ( .NOT. salsa_gases_from_chem )  THEN
828                DO  ig = 1, ngases_salsa
829                   CALL exchange_horiz( salsa_gas(ig)%conc_p, nbgp,                                &
830                                        alternative_communicator = communicator_salsa )
831                ENDDO
832             ENDIF
833          ENDIF
834
835          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' )
836
837!
838!--       Boundary conditions for the prognostic quantities (except of the
839!--       velocities at the outflow in case of a non-cyclic lateral wall) and
840!--       boundary conditions for module-specific variables
841          CALL module_interface_boundary_conditions
842
843!
844!--       Incrementing timestep counter
845          timestep_count = timestep_count + 1
846
847          CALL cpu_log( log_point(28), 'swap_timelevel', 'start' )
848!
849!--       Set the swap level for all modules
850          CALL module_interface_swap_timelevel( MOD( timestep_count, 2) )
851
852#if defined( __parallel )
853!
854!--       Set the swap level for steering the pmc data transfer
855          IF ( nested_run )  CALL pmci_set_swaplevel( MOD( timestep_count, 2) + 1 )  !> @todo: why the +1 ?
856#endif
857
858          CALL cpu_log( log_point(28), 'swap_timelevel', 'stop' )
859
860#if defined( __parallel )
861          IF ( nested_run )  THEN
862
863             CALL cpu_log( log_point(60), 'nesting', 'start' )
864!
865!--          Domain nesting. The data transfer subroutines pmci_parent_datatrans
866!--          and pmci_child_datatrans are called inside the wrapper
867!--          subroutine pmci_datatrans according to the control parameters
868!--          nesting_mode and nesting_datatransfer_mode.
869!--          TO_DO: why is nesting_mode given as a parameter here?
870             CALL pmci_datatrans( nesting_mode )
871
872             IF ( TRIM( nesting_mode ) == 'two-way' .OR.  nesting_mode == 'vertical' )  THEN
873
874                CALL cpu_log( log_point_s(92), 'exchange-horiz-nest', 'start' )
875!
876!--             Exchange_horiz is needed for all parent-domains after the
877!--             anterpolation
878                CALL exchange_horiz( u, nbgp )
879                CALL exchange_horiz( v, nbgp )
880                CALL exchange_horiz( w, nbgp )
881                IF ( .NOT. neutral )  CALL exchange_horiz( pt, nbgp )
882
883                IF ( humidity )  THEN
884
885                   CALL exchange_horiz( q, nbgp )
886
887                   IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
888                       CALL exchange_horiz( qc, nbgp )
889                       CALL exchange_horiz( nc, nbgp )
890                   ENDIF
891                   IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
892                       CALL exchange_horiz( qr, nbgp )
893                       CALL exchange_horiz( nr, nbgp )
894                   ENDIF
895                   IF ( bulk_cloud_model  .AND.  microphysics_ice_phase )  THEN
896                      CALL exchange_horiz( qi, nbgp )
897                      CALL exchange_horiz( ni, nbgp )
898                   ENDIF
899
900                ENDIF
901
902                IF ( passive_scalar )  CALL exchange_horiz( s, nbgp )
903
904                IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e, nbgp )
905
906                IF ( .NOT. constant_diffusion  .AND.  rans_mode  .AND.  rans_tke_e )  THEN
907                   CALL exchange_horiz( diss, nbgp )
908                ENDIF
909
910                IF ( air_chemistry )  THEN
911                   DO  n = 1, nvar
912                      CALL exchange_horiz( chem_species(n)%conc, nbgp,                             &
913                                           alternative_communicator = communicator_chem )
914                   ENDDO
915                ENDIF
916
917                IF ( salsa  .AND. time_since_reference_point >= skip_time_do_salsa )  THEN
918                   DO  ib = 1, nbins_aerosol
919                      CALL exchange_horiz( aerosol_number(ib)%conc, nbgp,                          &
920                                           alternative_communicator = communicator_salsa )
921                      DO  ic = 1, ncomponents_mass
922                         icc = ( ic - 1 ) * nbins_aerosol + ib
923                         CALL exchange_horiz( aerosol_mass(icc)%conc, nbgp,                        &
924                                              alternative_communicator = communicator_salsa )
925                      ENDDO
926                   ENDDO
927                   IF ( .NOT. salsa_gases_from_chem )  THEN
928                      DO  ig = 1, ngases_salsa
929                         CALL exchange_horiz( salsa_gas(ig)%conc, nbgp,                            &
930                                              alternative_communicator = communicator_salsa )
931                      ENDDO
932                   ENDIF
933                ENDIF
934                CALL cpu_log( log_point_s(92), 'exchange-horiz-nest', 'stop' )
935
936             ENDIF
937
938!
939!--          Set boundary conditions again after interpolation and anterpolation.
940             CALL pmci_boundary_conds
941
942             CALL cpu_log( log_point(60), 'nesting', 'stop' )
943
944          ENDIF
945#endif
946
947!
948!--       Temperature offset must be imposed at cyclic boundaries in x-direction
949!--       when a sloping surface is used
950          IF ( sloping_surface )  THEN
951             IF ( nxl ==  0 )  pt(:,:,nxlg:nxl-1) = pt(:,:,nxlg:nxl-1) - pt_slope_offset
952             IF ( nxr == nx )  pt(:,:,nxr+1:nxrg) = pt(:,:,nxr+1:nxrg) + pt_slope_offset
953          ENDIF
954
955!
956!--       Increase temperature pt(0) according to pt_surface_heating_rate (convert from K/h to K/s)
957          IF ( pt_surface_heating_rate /= 0.0_wp .AND. intermediate_timestep_count == 1 ) THEN
958             pt(0,:,:) = pt(0,:,:) + dt_3d * pt_surface_heating_rate / 3600.0_wp
959          ENDIF
960
961!
962!--       Impose a turbulent inflow using the recycling method
963          IF ( turbulent_inflow )  CALL inflow_turbulence
964
965!
966!--       Set values at outflow boundary using the special outflow condition
967          IF ( turbulent_outflow )  CALL outflow_turbulence
968
969!
970!--       Impose a random perturbation on the horizontal velocity field
971          IF ( create_disturbances  .AND.  ( call_psolver_at_all_substeps  .AND.                   &
972               intermediate_timestep_count == intermediate_timestep_count_max )                    &
973               .OR. ( .NOT. call_psolver_at_all_substeps  .AND.  intermediate_timestep_count == 1 ) ) &
974          THEN
975             time_disturb = time_disturb + dt_3d
976             IF ( time_disturb >= dt_disturb )  THEN
977                IF ( disturbance_energy_limit /= 0.0_wp  .AND.                                     &
978                     hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit )  THEN
979                   CALL disturb_field( 'u', tend, u )
980                   CALL disturb_field( 'v', tend, v )
981                ELSEIF ( ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )                                &
982                         .AND. .NOT. child_domain  .AND.  .NOT.  nesting_offline )                 &
983                THEN
984!
985!--                Runs with a non-cyclic lateral wall need perturbations
986!--                near the inflow throughout the whole simulation
987                   dist_range = 1
988                   CALL disturb_field( 'u', tend, u )
989                   CALL disturb_field( 'v', tend, v )
990                   dist_range = 0
991                ENDIF
992                time_disturb = time_disturb - dt_disturb
993             ENDIF
994          ENDIF
995
996!
997!--       Map forcing data derived from larger scale model onto domain
998!--       boundaries. Further, update geostrophic wind components.
999          IF ( nesting_offline  .AND.  intermediate_timestep_count ==                              &
1000                                       intermediate_timestep_count_max  )  THEN
1001!--          Determine interpolation factor before boundary conditions and geostrophic wind
1002!--          is updated.
1003             CALL nesting_offl_interpolation_factor
1004             CALL nesting_offl_bc
1005             CALL nesting_offl_geostrophic_wind
1006          ENDIF
1007!
1008!--       Impose a turbulent inflow using synthetic generated turbulence.
1009          IF ( use_syn_turb_gen  .AND.                                                             &
1010               intermediate_timestep_count == intermediate_timestep_count_max )  THEN
1011             CALL cpu_log( log_point(57), 'synthetic_turbulence_gen', 'start' )
1012             CALL stg_main
1013             CALL cpu_log( log_point(57), 'synthetic_turbulence_gen', 'stop' )
1014          ENDIF
1015!
1016!--       Ensure mass conservation. This need to be done after imposing
1017!--       synthetic turbulence and top boundary condition for pressure is set to
1018!--       Neumann conditions.
1019!--       Is this also required in case of Dirichlet?
1020          IF ( nesting_offline )  CALL nesting_offl_mass_conservation
1021!
1022!--       Reduce the velocity divergence via the equation for perturbation
1023!--       pressure.
1024          IF ( intermediate_timestep_count == 1  .OR.  call_psolver_at_all_substeps )  THEN
1025
1026#if defined( __parallel )
1027!
1028!--          Mass (volume) flux correction to ensure global mass conservation for child domains.
1029             IF ( child_domain )  THEN
1030                IF ( nesting_mode == 'vertical' )  THEN
1031                   CALL pmci_ensure_nest_mass_conservation_vertical
1032                ELSE
1033                   CALL pmci_ensure_nest_mass_conservation
1034                ENDIF
1035             ENDIF
1036#endif
1037             CALL pres
1038
1039          ENDIF
1040!
1041!--       Particle transport/physics with the Lagrangian particle model
1042!--       (only once during intermediate steps, because it uses an Euler-step)
1043!--       ### particle model should be moved before prognostic_equations, in order
1044!--       to regard droplet interactions directly
1045
1046          CALL module_interface_actions( 'after_pressure_solver' )
1047!
1048!--       Interaction of droplets with temperature and mixing ratio.
1049!--       Droplet condensation and evaporation is calculated within
1050!--       advec_particles.
1051!
1052!--       If required, compute liquid water content
1053          IF ( bulk_cloud_model )  THEN
1054             CALL calc_liquid_water_content
1055          ENDIF
1056!
1057!--       If required, compute virtual potential temperature
1058          IF ( humidity )  THEN
1059             CALL compute_vpt
1060          ENDIF
1061
1062!
1063!--       Compute the diffusion quantities
1064          IF ( .NOT. constant_diffusion )  THEN
1065
1066!
1067!--          Determine surface fluxes shf and qsws and surface values
1068!--          pt_surface and q_surface in dependence on data from external
1069!--          file LSF_DATA respectively
1070             IF ( ( large_scale_forcing .AND. lsf_surf ) .AND.                                     &
1071                 intermediate_timestep_count == intermediate_timestep_count_max )                  &
1072             THEN
1073                CALL ls_forcing_surf( simulated_time )
1074             ENDIF
1075
1076!
1077!--          First the vertical (and horizontal) fluxes in the surface
1078!--          (constant flux) layer are computed
1079             IF ( constant_flux_layer )  THEN
1080                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'start' )
1081                CALL surface_layer_fluxes
1082                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'stop' )
1083             ENDIF
1084!
1085!--          If required, solve the energy balance for the surface and run soil
1086!--          model. Call for horizontal as well as vertical surfaces
1087             IF ( land_surface .AND. time_since_reference_point >= skip_time_do_lsm)  THEN
1088
1089                CALL cpu_log( log_point(54), 'land_surface', 'start' )
1090!
1091!--             Call for horizontal upward-facing surfaces
1092                CALL lsm_energy_balance( .TRUE., -1 )
1093                CALL lsm_soil_model( .TRUE., -1, .TRUE. )
1094!
1095!--             Call for northward-facing surfaces
1096                CALL lsm_energy_balance( .FALSE., 0 )
1097                CALL lsm_soil_model( .FALSE., 0, .TRUE. )
1098!
1099!--             Call for southward-facing surfaces
1100                CALL lsm_energy_balance( .FALSE., 1 )
1101                CALL lsm_soil_model( .FALSE., 1, .TRUE. )
1102!
1103!--             Call for eastward-facing surfaces
1104                CALL lsm_energy_balance( .FALSE., 2 )
1105                CALL lsm_soil_model( .FALSE., 2, .TRUE. )
1106!
1107!--             Call for westward-facing surfaces
1108                CALL lsm_energy_balance( .FALSE., 3 )
1109                CALL lsm_soil_model( .FALSE., 3, .TRUE. )
1110
1111!
1112!--             At the end, set boundary conditons for potential temperature
1113!--             and humidity after running the land-surface model. This
1114!--             might be important for the nesting, where arrays are transfered.
1115                CALL lsm_boundary_condition
1116
1117
1118                CALL cpu_log( log_point(54), 'land_surface', 'stop' )
1119             ENDIF
1120!
1121!--          If required, solve the energy balance for urban surfaces and run
1122!--          the material heat model
1123             IF (urban_surface) THEN
1124                CALL cpu_log( log_point(74), 'urban_surface', 'start' )
1125
1126                CALL usm_surface_energy_balance( .FALSE. )
1127                IF ( usm_material_model )  THEN
1128                   CALL usm_green_heat_model
1129                   CALL usm_material_heat_model ( .FALSE. )
1130                ENDIF
1131
1132!
1133!--             At the end, set boundary conditons for potential temperature
1134!--             and humidity after running the urban-surface model. This
1135!--             might be important for the nesting, where arrays are transfered.
1136                CALL usm_boundary_condition
1137
1138                CALL cpu_log( log_point(74), 'urban_surface', 'stop' )
1139             ENDIF
1140!
1141!--          Compute the diffusion coefficients
1142             CALL cpu_log( log_point(17), 'diffusivities', 'start' )
1143             IF ( .NOT. humidity ) THEN
1144                IF ( ocean_mode )  THEN
1145                   CALL tcm_diffusivities( prho, prho_reference )
1146                ELSE
1147                   CALL tcm_diffusivities( pt, pt_reference )
1148                ENDIF
1149             ELSE
1150                CALL tcm_diffusivities( vpt, pt_reference )
1151             ENDIF
1152             CALL cpu_log( log_point(17), 'diffusivities', 'stop' )
1153
1154          ENDIF
1155
1156       ENDDO   ! Intermediate step loop
1157
1158!
1159!--    Will be used at some point by flow_statistics.
1160       !$ACC UPDATE &
1161       !$ACC HOST(sums_l_l(nzb:nzt+1,0:statistic_regions,0)) &
1162       !$ACC HOST(sums_us2_ws_l(nzb:nzt+1,0)) &
1163       !$ACC HOST(sums_wsus_ws_l(nzb:nzt+1,0)) &
1164       !$ACC HOST(sums_vs2_ws_l(nzb:nzt+1,0)) &
1165       !$ACC HOST(sums_wsvs_ws_l(nzb:nzt+1,0)) &
1166       !$ACC HOST(sums_ws2_ws_l(nzb:nzt+1,0)) &
1167       !$ACC HOST(sums_wspts_ws_l(nzb:nzt+1,0)) &
1168       !$ACC HOST(sums_wssas_ws_l(nzb:nzt+1,0)) &
1169       !$ACC HOST(sums_wsqs_ws_l(nzb:nzt+1,0)) &
1170       !$ACC HOST(sums_wsqcs_ws_l(nzb:nzt+1,0)) &
1171       !$ACC HOST(sums_wsqis_ws_l(nzb:nzt+1,0)) &
1172       !$ACC HOST(sums_wsqrs_ws_l(nzb:nzt+1,0)) &
1173       !$ACC HOST(sums_wsncs_ws_l(nzb:nzt+1,0)) &
1174       !$ACC HOST(sums_wsnis_ws_l(nzb:nzt+1,0)) &
1175       !$ACC HOST(sums_wsnrs_ws_l(nzb:nzt+1,0)) &
1176       !$ACC HOST(sums_wsss_ws_l(nzb:nzt+1,0)) &
1177       !$ACC HOST(sums_salsa_ws_l(nzb:nzt+1,0))
1178
1179!
1180!--    If required, calculate radiative fluxes and heating rates
1181       IF ( radiation  .AND.  time_since_reference_point > skip_time_do_radiation )  THEN
1182
1183          time_radiation = time_radiation + dt_3d
1184
1185          IF ( time_radiation >= dt_radiation  .OR.  force_radiation_call )  THEN
1186
1187             CALL cpu_log( log_point(50), 'radiation', 'start' )
1188
1189             IF ( .NOT. force_radiation_call )  THEN
1190                time_radiation = time_radiation - dt_radiation
1191             ENDIF
1192
1193!
1194!--          Adjust the current time to the time step of the radiation model.
1195!--          Needed since radiation is pre-calculated and stored only on apparent
1196!--          solar positions
1197             time_since_reference_point_save = time_since_reference_point
1198             time_since_reference_point = REAL( FLOOR( time_since_reference_point /             &
1199                                                       dt_radiation ), wp ) * dt_radiation
1200
1201             CALL radiation_control
1202
1203             IF ( ( urban_surface  .OR.  land_surface )  .AND.  radiation_interactions )  THEN
1204                CALL cpu_log( log_point_s(46), 'radiation_interaction', 'start' )
1205                CALL radiation_interaction
1206                CALL cpu_log( log_point_s(46), 'radiation_interaction', 'stop' )
1207             ENDIF
1208
1209!
1210!--          Return the current time to its original value
1211             time_since_reference_point = time_since_reference_point_save
1212
1213             CALL cpu_log( log_point(50), 'radiation', 'stop' )
1214
1215          ENDIF
1216       ENDIF
1217
1218
1219!
1220!-- 20200203 (ECC)
1221!-- allows for emission update mode in legacy mode as well as on-demand mode
1222!-- note that under on-demand mode emission update is no longer restricted to
1223!-- an hourly frequency, but whenever the simulation time corresponds to an
1224!-- inrement in emission timestamp value
1225
1226!
1227!-- If required, consider chemical emissions
1228
1229       IF  ( air_chemistry .AND. emissions_anthropogenic )  THEN
1230
1231          IF  ( emiss_read_legacy_mode )  THEN
1232!
1233!-- get hourly index and updates emission data when the hour is passed
1234
1235             CALL get_date_time( time_since_reference_point, hour=hour )
1236
1237             IF  ( hour_call_emis /= hour )   THEN
1238
1239                CALL chem_emissions_setup( chem_emis_att, chem_emis, n_matched_vars )
1240                hour_call_emis = hour
1241
1242             ENDIF
1243
1244          ELSE
1245
1246             CALL chem_emissions_update_on_demand
1247
1248          ENDIF
1249
1250       ENDIF
1251
1252
1253!
1254!--    If required, consider aerosol emissions for the salsa model
1255       IF ( salsa )  THEN
1256!
1257!--       Call emission routine to update emissions if needed
1258          CALL salsa_emission_update
1259
1260       ENDIF
1261!
1262!--    If required, calculate indoor temperature, waste heat, heat flux
1263!--    through wall, etc.
1264!--    dt_indoor steers the frequency of the indoor model calculations.
1265!--    Note, at first timestep indoor model is called, in order to provide
1266!--    a waste heat flux.
1267       IF ( indoor_model )  THEN
1268
1269          time_indoor = time_indoor + dt_3d
1270
1271          IF ( time_indoor >= dt_indoor  .OR.  current_timestep_number == 0 )  THEN
1272
1273             time_indoor = time_indoor - dt_indoor
1274
1275             CALL cpu_log( log_point(76), 'indoor_model', 'start' )
1276             CALL im_main_heatcool
1277             CALL cpu_log( log_point(76), 'indoor_model', 'stop' )
1278
1279          ENDIF
1280       ENDIF
1281!
1282!--    Increase simulation time and output times
1283       nr_timesteps_this_run      = nr_timesteps_this_run + 1
1284       current_timestep_number    = current_timestep_number + 1
1285       simulated_time             = simulated_time   + dt_3d
1286       time_since_reference_point = simulated_time - coupling_start_time
1287       simulated_time_chr         = time_to_string( time_since_reference_point )
1288
1289
1290       IF ( time_since_reference_point >= skip_time_data_output_av )  THEN
1291          time_do_av         = time_do_av       + dt_3d
1292       ENDIF
1293       IF ( time_since_reference_point >= skip_time_do2d_xy )  THEN
1294          time_do2d_xy       = time_do2d_xy     + dt_3d
1295       ENDIF
1296       IF ( time_since_reference_point >= skip_time_do2d_xz )  THEN
1297          time_do2d_xz       = time_do2d_xz     + dt_3d
1298       ENDIF
1299       IF ( time_since_reference_point >= skip_time_do2d_yz )  THEN
1300          time_do2d_yz       = time_do2d_yz     + dt_3d
1301       ENDIF
1302       IF ( time_since_reference_point >= skip_time_do3d    )  THEN
1303          time_do3d          = time_do3d        + dt_3d
1304       ENDIF
1305       DO  mid = 1, masks
1306          IF ( time_since_reference_point >= skip_time_domask(mid) )  THEN
1307             time_domask(mid)= time_domask(mid) + dt_3d
1308          ENDIF
1309       ENDDO
1310       IF ( time_since_reference_point >= skip_time_dosp )  THEN
1311          time_dosp       = time_dosp        + dt_3d
1312       ENDIF
1313       time_dots          = time_dots        + dt_3d
1314       IF ( .NOT. first_call_lpm )  THEN
1315          time_dopts      = time_dopts       + dt_3d
1316       ENDIF
1317       IF ( time_since_reference_point >= skip_time_dopr )  THEN
1318          time_dopr       = time_dopr        + dt_3d
1319       ENDIF
1320       time_dopr_listing  = time_dopr_listing + dt_3d
1321       time_run_control   = time_run_control + dt_3d
1322!
1323!--    Increment time-counter for surface output
1324       IF ( surface_output )  THEN
1325          IF ( time_since_reference_point >= skip_time_dosurf )  THEN
1326             time_dosurf    = time_dosurf + dt_3d
1327          ENDIF
1328          IF ( time_since_reference_point >= skip_time_dosurf_av )  THEN
1329             time_dosurf_av = time_dosurf_av + dt_3d
1330          ENDIF
1331       ENDIF
1332!
1333!--    Increment time-counter for virtual measurements
1334       IF ( virtual_measurement  .AND.  vm_time_start <= time_since_reference_point )  THEN
1335          time_virtual_measurement = time_virtual_measurement + dt_3d
1336       ENDIF
1337
1338!
1339!--    Increment time-counter for wind turbine data output
1340       IF ( wind_turbine )  THEN
1341          time_wtm = time_wtm + dt_3d
1342       ENDIF
1343
1344!
1345!--    In case of synthetic turbulence generation and parametrized turbulence
1346!--    information, update the time counter and if required, adjust the
1347!--    STG to new atmospheric conditions.
1348       IF ( use_syn_turb_gen  )  THEN
1349          IF ( parametrize_inflow_turbulence )  THEN
1350             time_stg_adjust = time_stg_adjust + dt_3d
1351             IF ( time_stg_adjust >= dt_stg_adjust )  THEN
1352                CALL cpu_log( log_point(57), 'synthetic_turbulence_gen', 'start' )
1353                CALL stg_adjust
1354                CALL cpu_log( log_point(57), 'synthetic_turbulence_gen', 'stop' )
1355             ENDIF
1356          ENDIF
1357          time_stg_call = time_stg_call + dt_3d
1358       ENDIF
1359
1360!
1361!--    Data exchange between coupled models
1362       IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
1363          time_coupling = time_coupling + dt_3d
1364
1365!
1366!--       In case of model termination initiated by the local model
1367!--       (terminate_coupled > 0), the coupler must be skipped because it would
1368!--       cause an MPI intercomminucation hang.
1369!--       If necessary, the coupler will be called at the beginning of the
1370!--       next restart run.
1371          DO WHILE ( time_coupling >= dt_coupling  .AND.  terminate_coupled == 0 )
1372             CALL surface_coupler
1373             time_coupling = time_coupling - dt_coupling
1374          ENDDO
1375       ENDIF
1376
1377!
1378!--    Biometeorology calculation of stationary thermal indices
1379!--    Todo (kanani): biometeorology needs own time_... treatment.
1380!--                   It might be that time_do2d_xy differs from time_do3d,
1381!--                   and then we might get trouble with the biomet output,
1382!--                   because we can have 2d and/or 3d biomet output!!
1383       IF ( biometeorology                                                                         &
1384            .AND. ( ( time_do3d >= dt_do3d  .AND.  time_since_reference_point >= skip_time_do3d )  &
1385                  .OR.                                                                             &
1386            ( time_do2d_xy >= dt_do2d_xy  .AND.  time_since_reference_point >= skip_time_do2d_xy ) &
1387                    ) )  THEN
1388!
1389!--       If required, do thermal comfort calculations
1390          IF ( thermal_comfort )  THEN
1391             CALL bio_calculate_thermal_index_maps ( .FALSE. )
1392          ENDIF
1393!
1394!--       If required, do UV exposure calculations
1395          IF ( uv_exposure )  THEN
1396             CALL bio_calculate_uv_exposure
1397          ENDIF
1398       ENDIF
1399
1400!
1401!--    Execute alle other module actions routunes
1402       CALL module_interface_actions( 'after_integration' )
1403
1404!
1405!--    If Galilei transformation is used, determine the distance that the
1406!--    model has moved so far
1407       IF ( galilei_transformation )  THEN
1408          advected_distance_x = advected_distance_x + u_gtrans * dt_3d
1409          advected_distance_y = advected_distance_y + v_gtrans * dt_3d
1410       ENDIF
1411
1412!
1413!--    Check, if restart is necessary (because cpu-time is expiring or
1414!--    because it is forced by user) and set stop flag
1415!--    This call is skipped if the remote model has already initiated a restart.
1416       IF ( .NOT. terminate_run )  CALL check_for_restart
1417
1418!
1419!--    Carry out statistical analysis and output at the requested output times.
1420!--    The MOD function is used for calculating the output time counters (like
1421!--    time_dopr) in order to regard a possible decrease of the output time
1422!--    interval in case of restart runs
1423
1424!
1425!--    Set a flag indicating that so far no statistics have been created
1426!--    for this time step
1427       flow_statistics_called = .FALSE.
1428
1429!
1430!--    If required, call flow_statistics for averaging in time
1431       IF ( averaging_interval_pr /= 0.0_wp  .AND.                                                 &
1432            ( dt_dopr - time_dopr ) <= averaging_interval_pr  .AND.                                &
1433            time_since_reference_point >= skip_time_dopr )  THEN
1434          time_dopr_av = time_dopr_av + dt_3d
1435          IF ( time_dopr_av >= dt_averaging_input_pr )  THEN
1436             do_sum = .TRUE.
1437             time_dopr_av = MOD( time_dopr_av, MAX( dt_averaging_input_pr, dt_3d ) )
1438          ENDIF
1439       ENDIF
1440       IF ( do_sum )  CALL flow_statistics
1441
1442!
1443!--    Sum-up 3d-arrays for later output of time-averaged 2d/3d/masked data
1444       IF ( averaging_interval /= 0.0_wp  .AND.                                                    &
1445            ( dt_data_output_av - time_do_av ) <= averaging_interval  .AND.                        &
1446            time_since_reference_point >= skip_time_data_output_av )                               &
1447       THEN
1448          time_do_sla = time_do_sla + dt_3d
1449          IF ( time_do_sla >= dt_averaging_input )  THEN
1450             IF ( current_timestep_number > timestep_number_at_prev_calc )                         &
1451                CALL doq_calculate
1452
1453             CALL sum_up_3d_data
1454             average_count_3d = average_count_3d + 1
1455             time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) )
1456          ENDIF
1457       ENDIF
1458!
1459!--    Average surface data
1460       IF ( surface_output )  THEN
1461          IF ( averaging_interval_surf /= 0.0_wp                                                   &
1462                .AND.  ( dt_dosurf_av - time_dosurf_av ) <= averaging_interval_surf                &
1463                .AND.  time_since_reference_point >= skip_time_dosurf_av )  THEN
1464             IF ( time_dosurf_av >= dt_averaging_input )  THEN
1465                CALL surface_data_output_averaging
1466                average_count_surf = average_count_surf + 1
1467             ENDIF
1468          ENDIF
1469       ENDIF
1470
1471!
1472!--    Calculate spectra for time averaging
1473       IF ( averaging_interval_sp /= 0.0_wp  .AND. ( dt_dosp - time_dosp ) <= averaging_interval_sp&
1474            .AND.  time_since_reference_point >= skip_time_dosp )  THEN
1475          time_dosp_av = time_dosp_av + dt_3d
1476          IF ( time_dosp_av >= dt_averaging_input_pr )  THEN
1477             CALL calc_spectra
1478             time_dosp_av = MOD( time_dosp_av, MAX( dt_averaging_input_pr, dt_3d ) )
1479          ENDIF
1480       ENDIF
1481
1482!
1483!--    Call flight module and output data
1484       IF ( virtual_flight )  THEN
1485          CALL flight_measurement
1486          CALL data_output_flight
1487       ENDIF
1488!
1489!--    Take virtual measurements
1490       IF ( virtual_measurement  .AND.  time_virtual_measurement >= dt_virtual_measurement         &
1491                                 .AND.  vm_time_start <= time_since_reference_point )  THEN
1492          CALL vm_sampling
1493          CALL vm_data_output
1494          time_virtual_measurement = MOD(      time_virtual_measurement,                           &
1495                                          MAX( dt_virtual_measurement, dt_3d ) )
1496       ENDIF
1497
1498!
1499!--    Output wind turbine data
1500       IF ( wind_turbine  .AND.  time_wtm >= dt_data_output_wtm )  THEN
1501          CALL wtm_data_output
1502          time_wtm = MOD( time_wtm, MAX( dt_data_output_wtm, dt_3d ) )
1503       ENDIF
1504
1505!
1506!--    Profile output (ASCII) on file
1507       IF ( time_dopr_listing >= dt_dopr_listing )  THEN
1508          CALL print_1d
1509          time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, dt_3d ) )
1510       ENDIF
1511
1512!
1513!--    Graphic output for PROFIL
1514       IF ( time_dopr >= dt_dopr  .AND.  time_since_reference_point >= skip_time_dopr )  THEN
1515          IF ( dopr_n /= 0 )  CALL data_output_profiles
1516          time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) )
1517          time_dopr_av = 0.0_wp    ! due to averaging (see above)
1518       ENDIF
1519
1520!
1521!--    Graphic output for time series
1522       IF ( time_dots >= dt_dots )  THEN
1523          CALL data_output_tseries
1524          time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) )
1525       ENDIF
1526
1527!
1528!--    Output of spectra (formatted for use with PROFIL), in case of no
1529!--    time averaging, spectra has to be calculated before
1530       IF ( time_dosp >= dt_dosp  .AND.  time_since_reference_point >= skip_time_dosp )  THEN
1531          IF ( average_count_sp == 0 )  CALL calc_spectra
1532          CALL data_output_spectra
1533          time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) )
1534       ENDIF
1535
1536!
1537!--    2d-data output (cross-sections)
1538       IF ( time_do2d_xy >= dt_do2d_xy  .AND.  time_since_reference_point >= skip_time_do2d_xy )  THEN
1539          IF ( current_timestep_number > timestep_number_at_prev_calc )                            &
1540             CALL doq_calculate
1541
1542          CALL data_output_2d( 'xy', 0 )
1543          time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) )
1544       ENDIF
1545       IF ( time_do2d_xz >= dt_do2d_xz  .AND.  time_since_reference_point >= skip_time_do2d_xz )  THEN
1546          IF ( current_timestep_number > timestep_number_at_prev_calc )                            &
1547
1548             CALL doq_calculate
1549          CALL data_output_2d( 'xz', 0 )
1550          time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) )
1551       ENDIF
1552       IF ( time_do2d_yz >= dt_do2d_yz  .AND.  time_since_reference_point >= skip_time_do2d_yz )  THEN
1553          IF ( current_timestep_number > timestep_number_at_prev_calc )                            &
1554             CALL doq_calculate
1555
1556          CALL data_output_2d( 'yz', 0 )
1557          time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) )
1558       ENDIF
1559
1560!
1561!--    3d-data output (volume data)
1562       IF ( time_do3d >= dt_do3d  .AND.  time_since_reference_point >= skip_time_do3d )  THEN
1563          IF ( current_timestep_number > timestep_number_at_prev_calc )                            &
1564             CALL doq_calculate
1565
1566          CALL data_output_3d( 0 )
1567          time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) )
1568       ENDIF
1569
1570!
1571!--    Masked data output
1572       DO  mid = 1, masks
1573          IF ( time_domask(mid) >= dt_domask(mid)                                                  &
1574               .AND.  time_since_reference_point >= skip_time_domask(mid) )  THEN
1575             IF ( current_timestep_number > timestep_number_at_prev_calc )                         &
1576                CALL doq_calculate
1577
1578             CALL data_output_mask( 0, mid )
1579             time_domask(mid) = MOD( time_domask(mid), MAX( dt_domask(mid), dt_3d ) )
1580          ENDIF
1581       ENDDO
1582
1583!
1584!--    Output of time-averaged 2d/3d/masked data
1585       IF ( time_do_av >= dt_data_output_av                                                        &
1586            .AND.  time_since_reference_point >= skip_time_data_output_av )  THEN
1587          CALL average_3d_data
1588!
1589!--       Udate thermal comfort indices based on updated averaged input
1590          IF ( biometeorology  .AND.  thermal_comfort )  THEN
1591             CALL bio_calculate_thermal_index_maps ( .TRUE. )
1592          ENDIF
1593          CALL data_output_2d( 'xy', 1 )
1594          CALL data_output_2d( 'xz', 1 )
1595          CALL data_output_2d( 'yz', 1 )
1596          CALL data_output_3d( 1 )
1597          DO  mid = 1, masks
1598             CALL data_output_mask( 1, mid )
1599          ENDDO
1600          time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) )
1601       ENDIF
1602!
1603!--    Output of surface data, instantaneous and averaged data
1604       IF ( surface_output )  THEN
1605          IF ( time_dosurf >= dt_dosurf  .AND.  time_since_reference_point >= skip_time_dosurf )  THEN
1606             CALL surface_data_output( 0 )
1607             time_dosurf = MOD( time_dosurf, MAX( dt_dosurf, dt_3d ) )
1608          ENDIF
1609          IF ( time_dosurf_av >= dt_dosurf_av  .AND.  time_since_reference_point >= skip_time_dosurf_av )  THEN
1610             CALL surface_data_output( 1 )
1611             time_dosurf_av = MOD( time_dosurf_av, MAX( dt_dosurf_av, dt_3d ) )
1612          ENDIF
1613       ENDIF
1614
1615!
1616!--    Output of particle time series
1617       IF ( particle_advection )  THEN
1618          IF ( time_dopts >= dt_dopts  .OR.                                                        &
1619               ( time_since_reference_point >= particle_advection_start  .AND.                     &
1620                 first_call_lpm ) )  THEN
1621             CALL lpm_data_output_ptseries
1622             time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
1623          ENDIF
1624       ENDIF
1625
1626!
1627!--    If required, set the heat flux for the next time step to a random value
1628       IF ( constant_heatflux  .AND.  random_heatflux )  THEN
1629          IF ( surf_def_h(0)%ns >= 1 )  THEN
1630             CALL cpu_log( log_point(23), 'disturb_heatflux', 'start' )
1631             CALL disturb_heatflux( surf_def_h(0) )
1632             CALL cpu_log( log_point(23), 'disturb_heatflux', 'stop' )
1633          ENDIF
1634          IF ( surf_lsm_h%ns    >= 1 )  THEN
1635             CALL cpu_log( log_point(23), 'disturb_heatflux', 'start' )
1636             CALL disturb_heatflux( surf_lsm_h    )
1637             CALL cpu_log( log_point(23), 'disturb_heatflux', 'stop' )
1638          ENDIF
1639          IF ( surf_usm_h%ns    >= 1 )  THEN
1640             CALL cpu_log( log_point(23), 'disturb_heatflux', 'start' )
1641             CALL disturb_heatflux( surf_usm_h    )
1642             CALL cpu_log( log_point(23), 'disturb_heatflux', 'stop' )
1643          ENDIF
1644       ENDIF
1645
1646!
1647!--    Execute alle other module actions routunes
1648       CALL module_interface_actions( 'after_timestep' )
1649
1650!
1651!--    Determine size of next time step. Save timestep dt_3d because it is
1652!--    newly calculated in routine timestep, but required further below for
1653!--    steering the run control output interval
1654       dt_3d_old = dt_3d
1655       CALL timestep
1656
1657#if defined( __parallel )
1658!
1659!--    Synchronize the timestep in case of nested run.
1660       IF ( nested_run )  THEN
1661!
1662!--       Synchronize by unifying the time step.
1663!--       Global minimum of all time-steps is used for all.
1664          CALL pmci_synchronize
1665       ENDIF
1666#endif
1667
1668!
1669!--    Computation and output of run control parameters.
1670!--    This is also done whenever perturbations have been imposed
1671       IF ( time_run_control >= dt_run_control  .OR.                                               &
1672            timestep_scheme(1:5) /= 'runge'  .OR.  disturbance_created )                           &
1673       THEN
1674          CALL run_control
1675          IF ( time_run_control >= dt_run_control )  THEN
1676             time_run_control = MOD( time_run_control, MAX( dt_run_control, dt_3d_old ) )
1677          ENDIF
1678       ENDIF
1679
1680!
1681!--    Output elapsed simulated time in form of a progress bar on stdout
1682       IF ( myid == 0 )  CALL output_progress_bar
1683
1684       CALL cpu_log( log_point_s(10), 'timesteps', 'stop' )
1685
1686
1687    ENDDO   ! time loop
1688
1689#if defined( _OPENACC )
1690    CALL exit_surface_arrays
1691#endif
1692!$ACC END DATA
1693!$ACC END DATA
1694!$ACC END DATA
1695!$ACC END DATA
1696!$ACC END DATA
1697!$ACC END DATA
1698!$ACC END DATA
1699!$ACC END DATA
1700
1701    IF ( myid == 0 )  CALL finish_progress_bar
1702
1703    CALL location_message( 'atmosphere (and/or ocean) time-stepping', 'finished' )
1704
1705 END SUBROUTINE time_integration
Note: See TracBrowser for help on using the repository browser.