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

Last change on this file since 4462 was 4457, checked in by raasch, 5 years ago

ghost point exchange modularized, bugfix for wrong 2d-exchange

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