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

Last change on this file since 4540 was 4521, checked in by schwenkel, 4 years ago

add error number

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