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

Last change on this file since 3182 was 3182, checked in by suehring, 6 years ago

New Inifor features: grid stretching, improved command-interface, support start dates in different formats in both YYYYMMDD and YYYYMMDDHH, Ability to manually control input file prefixes (--radiation-prefix, --soil-preifx, --flow-prefix, --soilmoisture-prefix) for compatiblity with DWD forcast naming scheme; GNU-style short and long option; Prepared output of large-scale forcing profiles (no computation yet); Added preprocessor flag netcdf4 to switch output format between netCDF 3 and 4; Updated netCDF variable names and attributes to comply with PIDS v1.9; Inifor bugfixes: Improved compatibility with older Intel Intel compilers by avoiding implicit array allocation; Added origin_lon/_lat values and correct reference time in dynamic driver global attributes; corresponding PALM changes: adjustments to revised Inifor; variables names in dynamic driver adjusted; enable geostrophic forcing also in offline nested mode; variable names in LES-LES and COSMO offline nesting changed; lateral boundary flags for nesting, in- and outflow conditions renamed

  • Property svn:keywords set to Id
File size: 53.6 KB
Line 
1!> @file time_integration.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2018 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22! Replace simulated_time by time_since_reference_point in COSMO nesting mode.
23! Rename subroutines and variables in COSMO nesting mode
24!
25! Former revisions:
26! -----------------
27! $Id: time_integration.f90 3182 2018-07-27 13:36:03Z suehring $
28! Added multi agent system
29!
30! 3042 2018-05-25 10:44:37Z schwenkel
31! Changed the name specific humidity to mixing ratio
32!
33! 3040 2018-05-25 10:22:08Z schwenkel
34! Fixed bug in IF statement
35! Ensure that the time when calling the radiation to be the time step of the
36! pre-calculated time when first calculate the positions of the sun
37!
38! 3004 2018-04-27 12:33:25Z Giersch
39! First call of flow_statistics has been removed. It is already called in
40! run_control itself
41!
42! 2984 2018-04-18 11:51:30Z hellstea
43! CALL pmci_ensure_nest_mass_conservation is removed (so far only commented out)
44! as seemingly unnecessary.
45!
46! 2941 2018-04-03 11:54:58Z kanani
47! Deduct spinup_time from RUN_CONTROL output of main 3d run
48! (use time_since_reference_point instead of simulated_time)
49!
50! 2938 2018-03-27 15:52:42Z suehring
51! Nesting of dissipation rate in case of RANS mode and TKE-e closure is applied
52!
53! 2936 2018-03-27 14:49:27Z suehring
54! Little formatting adjustment.
55!
56! 2817 2018-02-19 16:32:21Z knoop
57! Preliminary gust module interface implemented
58!
59! 2801 2018-02-14 16:01:55Z thiele
60! Changed lpm from subroutine to module.
61! Introduce particle transfer in nested models.
62!
63! 2776 2018-01-31 10:44:42Z Giersch
64! Variable use_synthetic_turbulence_generator has been abbreviated
65!
66! 2773 2018-01-30 14:12:54Z suehring
67! - Nesting for chemical species
68!
69! 2766 2018-01-22 17:17:47Z kanani
70! Removed preprocessor directive __chem
71!
72! 2718 2018-01-02 08:49:38Z maronga
73! Corrected "Former revisions" section
74!
75! 2696 2017-12-14 17:12:51Z kanani
76! - Change in file header (GPL part)
77! - Implementation of uv exposure model (FK)
78! - Moved vnest_boundary_conds_khkm from tcm_diffusivities to here (TG)
79! - renamed diffusivities to tcm_diffusivities (TG)
80! - implement prognostic equation for diss (TG)
81! - Moved/commented CALL to chem_emissions (FK)
82! - Added CALL to chem_emissions (FK)
83! - Implementation of chemistry module (FK)
84! - Calls for setting boundary conditions in USM and LSM (MS)
85! - Large-scale forcing with larger-scale models implemented (MS)
86! - Rename usm_radiation into radiation_interactions; merge with branch
87!   radiation (MS)
88! - added call for usm_green_heat_model for green building surfaces (RvT)
89! - added call for usm_temperature_near_surface for use in indoor model (RvT)
90!
91! 2617 2017-11-16 12:47:24Z suehring
92! Bugfix, assure that the reference state does not become zero.
93!
94! 2563 2017-10-19 15:36:10Z Giersch
95! Variable wind_turbine moved to module control_parameters
96!
97! 2365 2017-08-21 14:59:59Z kanani
98! Vertical grid nesting implemented (SadiqHuq)
99!
100! 2320 2017-07-21 12:47:43Z suehring
101! Set bottom boundary conditions after nesting interpolation and anterpolation
102!
103! 2299 2017-06-29 10:14:38Z maronga
104! Call of soil model adjusted
105!
106! 2292 2017-06-20 09:51:42Z schwenkel
107! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
108! includes two more prognostic equations for cloud drop concentration (nc) 
109! and cloud water content (qc).
110!
111! 2271 2017-06-09 12:34:55Z sward
112! Start timestep message changed
113!
114! 2259 2017-06-08 09:09:11Z gronemeier
115! Implemented synthetic turbulence generator
116!
117! 2233 2017-05-30 18:08:54Z suehring
118!
119! 2232 2017-05-30 17:47:52Z suehring
120! Adjustments to new topography and surface concept
121! Modify passed parameters for disturb_field
122!
123! 2178 2017-03-17 11:07:39Z hellstea
124! Setting perturbations at all times near inflow boundary is removed
125! in case of nested boundaries
126!
127! 2174 2017-03-13 08:18:57Z maronga
128! Added support for nesting with cloud microphysics
129!
130! 2118 2017-01-17 16:38:49Z raasch
131! OpenACC directives and related code removed
132!
133! 2050 2016-11-08 15:00:55Z gronemeier
134! Implement turbulent outflow condition
135!
136! 2031 2016-10-21 15:11:58Z knoop
137! renamed variable rho to rho_ocean
138!
139! 2011 2016-09-19 17:29:57Z kanani
140! Flag urban_surface is now defined in module control_parameters,
141! removed commented CALLs of global_min_max.
142!
143! 2007 2016-08-24 15:47:17Z kanani
144! Added CALLs for new urban surface model
145!
146! 2000 2016-08-20 18:09:15Z knoop
147! Forced header and separation lines into 80 columns
148!
149! 1976 2016-07-27 13:28:04Z maronga
150! Simplified calls to radiation model
151!
152! 1960 2016-07-12 16:34:24Z suehring
153! Separate humidity and passive scalar
154!
155! 1957 2016-07-07 10:43:48Z suehring
156! flight module added
157!
158! 1919 2016-05-27 14:51:23Z raasch
159! Initial version of purely vertical nesting introduced.
160!
161! 1918 2016-05-27 14:35:57Z raasch
162! determination of time step moved to the end of the time step loop,
163! the first time step is now always calculated before the time step loop (i.e.
164! also in case of restart runs)
165!
166! 1914 2016-05-26 14:44:07Z witha
167! Added call for wind turbine model
168!
169! 1878 2016-04-19 12:30:36Z hellstea
170! Synchronization for nested runs rewritten
171!
172! 1853 2016-04-11 09:00:35Z maronga
173! Adjusted for use with radiation_scheme = constant
174!
175! 1849 2016-04-08 11:33:18Z hoffmann
176! Adapted for modularization of microphysics
177!
178! 1833 2016-04-07 14:23:03Z raasch
179! spectrum renamed spectra_mod, spectra related variables moved to spectra_mod
180!
181! 1831 2016-04-07 13:15:51Z hoffmann
182! turbulence renamed collision_turbulence
183!
184! 1822 2016-04-07 07:49:42Z hoffmann
185! icloud_scheme replaced by microphysics_*
186!
187! 1808 2016-04-05 19:44:00Z raasch
188! output message in case unscheduled radiation calls removed
189!
190! 1797 2016-03-21 16:50:28Z raasch
191! introduction of different datatransfer modes
192!
193! 1791 2016-03-11 10:41:25Z raasch
194! call of pmci_update_new removed
195!
196! 1786 2016-03-08 05:49:27Z raasch
197! +module spectrum
198!
199! 1783 2016-03-06 18:36:17Z raasch
200! switch back of netcdf data format for mask output moved to the mask output
201! routine
202!
203! 1781 2016-03-03 15:12:23Z raasch
204! some pmc calls removed at the beginning (before timeloop),
205! pmc initialization moved to the main program
206!
207! 1764 2016-02-28 12:45:19Z raasch
208! PMC_ACTIVE flags removed,
209! bugfix: nest synchronization after first call of timestep
210!
211! 1762 2016-02-25 12:31:13Z hellstea
212! Introduction of nested domain feature
213!
214! 1736 2015-12-04 08:56:33Z raasch
215! no perturbations added to total domain if energy limit has been set zero
216!
217! 1691 2015-10-26 16:17:44Z maronga
218! Added option for spin-ups without land surface and radiation models. Moved calls
219! for radiation and lan surface schemes.
220!
221! 1682 2015-10-07 23:56:08Z knoop
222! Code annotations made doxygen readable
223!
224! 1671 2015-09-25 03:29:37Z raasch
225! bugfix: ghostpoint exchange for array diss in case that sgs velocities are used
226! for particles
227!
228! 1585 2015-04-30 07:05:52Z maronga
229! Moved call of radiation scheme. Added support for RRTM
230!
231! 1551 2015-03-03 14:18:16Z maronga
232! Added interface for different radiation schemes.
233!
234! 1496 2014-12-02 17:25:50Z maronga
235! Added calls for the land surface model and radiation scheme
236!
237! 1402 2014-05-09 14:25:13Z raasch
238! location messages modified
239!
240! 1384 2014-05-02 14:31:06Z raasch
241! location messages added
242!
243! 1380 2014-04-28 12:40:45Z heinze
244! CALL of nudge_ref added
245! bc_pt_t_val and bc_q_t_val are updated in case nudging is used
246!
247! 1365 2014-04-22 15:03:56Z boeske
248! Reset sums_ls_l to zero at each timestep
249! +sums_ls_l
250! Calculation of reference state (previously in subroutine calc_mean_profile)
251
252! 1342 2014-03-26 17:04:47Z kanani
253! REAL constants defined as wp-kind
254!
255! 1320 2014-03-20 08:40:49Z raasch
256! ONLY-attribute added to USE-statements,
257! kind-parameters added to all INTEGER and REAL declaration statements,
258! kinds are defined in new module kinds,
259! old module precision_kind is removed,
260! revision history before 2012 removed,
261! comment fields (!:) to be used for variable explanations added to
262! all variable declaration statements
263! 1318 2014-03-17 13:35:16Z raasch
264! module interfaces removed
265!
266! 1308 2014-03-13 14:58:42Z fricke
267! +netcdf_data_format_save
268! For masked data, parallel netcdf output is not tested so far, hence
269! netcdf_data_format is switched back to non-paralell output.
270!
271! 1276 2014-01-15 13:40:41Z heinze
272! Use LSF_DATA also in case of Dirichlet bottom boundary condition for scalars
273!
274! 1257 2013-11-08 15:18:40Z raasch
275! acc-update-host directive for timestep removed
276!
277! 1241 2013-10-30 11:36:58Z heinze
278! Generalize calc_mean_profile for wider use
279! Determine shf and qsws in dependence on data from LSF_DATA
280! Determine ug and vg in dependence on data from LSF_DATA
281! 1221 2013-09-10 08:59:13Z raasch
282! host update of arrays before timestep is called
283!
284! 1179 2013-06-14 05:57:58Z raasch
285! mean profiles for reference state are only calculated if required,
286! small bugfix for background communication
287!
288! 1171 2013-05-30 11:27:45Z raasch
289! split of prognostic_equations deactivated (comment lines), for the time being
290!
291! 1128 2013-04-12 06:19:32Z raasch
292! asynchronous transfer of ghost point data realized for acc-optimized version:
293! prognostic_equations are first called two times for those points required for
294! the left-right and north-south exchange, respectively, and then for the
295! remaining points,
296! those parts requiring global communication moved from prognostic_equations to
297! here
298!
299! 1115 2013-03-26 18:16:16Z hoffmann
300! calculation of qr and nr is restricted to precipitation
301!
302! 1113 2013-03-10 02:48:14Z raasch
303! GPU-porting of boundary conditions,
304! openACC directives updated
305! formal parameter removed from routine boundary_conds
306!
307! 1111 2013-03-08 23:54:10Z raasch
308! +internal timestep counter for cpu statistics added,
309! openACC directives updated
310!
311! 1092 2013-02-02 11:24:22Z raasch
312! unused variables removed
313!
314! 1065 2012-11-22 17:42:36Z hoffmann
315! exchange of diss (dissipation rate) in case of turbulence = .TRUE. added
316!
317! 1053 2012-11-13 17:11:03Z hoffmann
318! exchange of ghost points for nr, qr added
319!
320! 1036 2012-10-22 13:43:42Z raasch
321! code put under GPL (PALM 3.9)
322!
323! 1019 2012-09-28 06:46:45Z raasch
324! non-optimized version of prognostic_equations removed
325!
326! 1015 2012-09-27 09:23:24Z raasch
327! +call of prognostic_equations_acc
328!
329! 1001 2012-09-13 14:08:46Z raasch
330! all actions concerning leapfrog- and upstream-spline-scheme removed
331!
332! 849 2012-03-15 10:35:09Z raasch
333! advec_particles renamed lpm, first_call_advec_particles renamed first_call_lpm
334!
335! 825 2012-02-19 03:03:44Z raasch
336! wang_collision_kernel renamed wang_kernel
337!
338! Revision 1.1  1997/08/11 06:19:04  raasch
339! Initial revision
340!
341!
342! Description:
343! ------------
344!> Integration in time of the model equations, statistical analysis and graphic
345!> output
346!------------------------------------------------------------------------------!
347 SUBROUTINE time_integration
348 
349
350    USE advec_ws,                                                              &
351        ONLY:  ws_statistics
352
353    USE arrays_3d,                                                             &
354        ONLY:  diss, diss_p, dzu, e, e_p, nc, nc_p, nr, nr_p, prho, pt, pt_p, pt_init, &
355               q_init, q, qc, qc_p, ql, ql_c, ql_v, ql_vp, qr, qr_p, q_p,      &
356               ref_state, rho_ocean, s, s_p, sa_p, tend, u, u_p, v, vpt,       &
357               v_p, w, w_p
358
359    USE calc_mean_profile_mod,                                                 &
360        ONLY:  calc_mean_profile
361
362    USE chemistry_model_mod,                                                   &
363        ONLY:  chem_emissions, chem_species
364
365    USE chem_modules,                                                          &
366        ONLY:  nspec 
367
368    USE control_parameters,                                                    &
369        ONLY:  advected_distance_x, advected_distance_y, air_chemistry,        &
370               average_count_3d, averaging_interval, averaging_interval_pr,    &
371               bc_lr_cyc, bc_ns_cyc, bc_pt_t_val, bc_q_t_val,                  &
372               call_psolver_at_all_substeps,  child_domain, cloud_droplets,    &
373               cloud_physics, constant_flux_layer, constant_heatflux,          &
374               create_disturbances, dopr_n, constant_diffusion, coupling_mode, &
375               coupling_start_time, current_timestep_number,                   &
376               disturbance_created, disturbance_energy_limit, dist_range,      &
377               do_sum, dt_3d, dt_averaging_input, dt_averaging_input_pr,       &
378               dt_coupling, dt_data_output_av, dt_disturb, dt_do2d_xy,         &
379               dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_domask,dt_dopts, dt_dopr,   &
380               dt_dopr_listing, dt_dots, dt_dvrp, dt_run_control, end_time,    &
381               first_call_lpm, first_call_mas, galilei_transformation,         &
382               humidity, intermediate_timestep_count,                          &
383               intermediate_timestep_count_max,                                &
384               land_surface, large_scale_forcing,                              &
385               loop_optimization, lsf_surf, lsf_vert, masks,                   &
386               microphysics_morrison, microphysics_seifert, mid,               &
387               nesting_offline, neutral, nr_timesteps_this_run, nudging,       &
388               ocean, passive_scalar, prho_reference, pt_reference,            &
389               pt_slope_offset, random_heatflux, rans_mode,                    &
390               rans_tke_e, run_coupled, simulated_time, simulated_time_chr,    &
391               skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,        &
392               skip_time_do3d, skip_time_domask, skip_time_dopr,               &
393               skip_time_data_output_av, sloping_surface,                      &
394               stop_dt, terminate_coupled, terminate_run, timestep_scheme,     &
395               time_coupling, time_do2d_xy, time_do2d_xz, time_do2d_yz,        &
396               time_do3d, time_domask, time_dopr, time_dopr_av,                &
397               time_dopr_listing, time_dopts, time_dosp, time_dosp_av,         &
398               time_dots, time_do_av, time_do_sla, time_disturb, time_dvrp,    &
399               time_run_control, time_since_reference_point,                   &
400               turbulent_inflow, turbulent_outflow, urban_surface,             &
401               use_initial_profile_as_reference,                               &
402               use_single_reference_value, uv_exposure, u_gtrans, v_gtrans,    &
403               virtual_flight, wind_turbine, ws_scheme_mom, ws_scheme_sca
404
405    USE cpulog,                                                                &
406        ONLY:  cpu_log, log_point, log_point_s
407
408    USE flight_mod,                                                            &
409        ONLY:  flight_measurement
410
411    USE gust_mod,                                                              &
412        ONLY:  gust_actions, gust_module_enabled
413
414    USE indices,                                                               &
415        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt
416
417    USE interaction_droplets_ptq_mod,                                          &
418        ONLY:  interaction_droplets_ptq
419
420    USE interfaces
421
422    USE kinds
423
424    USE land_surface_model_mod,                                                &
425        ONLY:  lsm_boundary_condition, lsm_energy_balance, lsm_soil_model,     &
426               skip_time_do_lsm
427
428    USE lsf_nudging_mod,                                                       &
429        ONLY:  calc_tnudge, ls_forcing_surf, ls_forcing_vert, nudge_ref,       &
430               lsf_nesting_offline, lsf_nesting_offline_mass_conservation
431
432    USE microphysics_mod,                                                      &
433        ONLY: collision_turbulence
434
435    USE netcdf_data_input_mod,                                                 &
436        ONLY:  nest_offl, netcdf_data_input_lsf
437
438    USE multi_agent_system_mod,                                                &
439        ONLY:  agents_active, multi_agent_system, multi_agent_system_start
440
441    USE particle_attributes,                                                   &
442        ONLY:  particle_advection, particle_advection_start,                   &
443               use_sgs_for_particles, wang_kernel
444
445    USE pegrid
446
447    USE pmc_interface,                                                         &
448        ONLY:  nested_run, nesting_mode, pmci_boundary_conds, pmci_datatrans,  &
449               pmci_ensure_nest_mass_conservation, pmci_synchronize
450
451    USE progress_bar,                                                          &
452        ONLY:  finish_progress_bar, output_progress_bar
453
454    USE prognostic_equations_mod,                                              &
455        ONLY:  prognostic_equations_cache, prognostic_equations_vector
456
457    USE radiation_model_mod,                                                   &
458        ONLY: dt_radiation, force_radiation_call, radiation, radiation_control,&
459              radiation_interaction, radiation_interactions,                   &
460              skip_time_do_radiation, time_radiation
461
462    USE spectra_mod,                                                           &
463        ONLY: average_count_sp, averaging_interval_sp, calc_spectra, dt_dosp,  &
464              skip_time_dosp
465
466    USE statistics,                                                            &
467        ONLY:  flow_statistics_called, hom, pr_palm, sums_ls_l, u_max,         &
468               u_max_ijk, v_max, v_max_ijk, w_max, w_max_ijk
469
470    USE surface_layer_fluxes_mod,                                              &
471        ONLY:  surface_layer_fluxes
472
473    USE surface_mod,                                                           &
474        ONLY:  surf_def_h, surf_lsm_h, surf_usm_h
475
476    USE turbulence_closure_mod,                                                &
477        ONLY:  tcm_diffusivities, production_e_init
478
479    USE urban_surface_mod,                                                     &
480        ONLY:  usm_boundary_condition, usm_material_heat_model,                &
481               usm_material_model,                                             &
482               usm_surface_energy_balance, usm_green_heat_model,               &
483               usm_temperature_near_surface
484
485    USE synthetic_turbulence_generator_mod,                                    &
486        ONLY:  stg_main, use_syn_turb_gen
487
488    USE user_actions_mod,                                                      &
489        ONLY:  user_actions
490
491    USE uv_exposure_model_mod,                                                 &
492        ONLY:  uvem_calc_exposure
493
494    USE wind_turbine_model_mod,                                                &
495        ONLY:  wtm_forces
496
497    USE lpm_mod,                                                               &
498        ONLY:  lpm
499
500    USE vertical_nesting_mod,                                                  &
501        ONLY:  vnested, vnest_anterpolate, vnest_anterpolate_e,                &
502               vnest_boundary_conds, vnest_boundary_conds_khkm,                & 
503               vnest_deallocate, vnest_init, vnest_init_fine,                  &
504               vnest_start_time
505
506    IMPLICIT NONE
507
508    CHARACTER (LEN=9) ::  time_to_string          !<
509    INTEGER(iwp)      ::  it
510    INTEGER(iwp)      ::  lsp
511    INTEGER(iwp)      ::  n
512
513
514    REAL(wp) ::  dt_3d_old  !< temporary storage of timestep to be used for
515                            !< steering of run control output interval
516    REAL(wp) ::  tsrp_org   !< original value of time_since_reference_point
517!
518!-- At beginning determine the first time step
519    CALL timestep
520!
521!-- Synchronize the timestep in case of nested run.
522    IF ( nested_run )  THEN
523!
524!--    Synchronization by unifying the time step.
525!--    Global minimum of all time-steps is used for all.
526       CALL pmci_synchronize
527    ENDIF
528
529!
530!-- Determine and print out the run control quantities before the first time
531!-- step of this run. For the initial run, some statistics (e.g. divergence)
532!-- need to be determined first --> CALL flow_statistics at the beginning of
533!-- run_control
534    CALL run_control
535!
536!-- Data exchange between coupled models in case that a call has been omitted
537!-- at the end of the previous run of a job chain.
538    IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled .AND. .NOT. vnested)  THEN
539!
540!--    In case of model termination initiated by the local model the coupler
541!--    must not be called because this would again cause an MPI hang.
542       DO WHILE ( time_coupling >= dt_coupling  .AND.  terminate_coupled == 0 )
543          CALL surface_coupler
544          time_coupling = time_coupling - dt_coupling
545       ENDDO
546       IF (time_coupling == 0.0_wp  .AND.                                      &
547           time_since_reference_point < dt_coupling )                          &
548       THEN
549          time_coupling = time_since_reference_point
550       ENDIF
551    ENDIF
552
553#if defined( __dvrp_graphics )
554!
555!-- Time measurement with dvrp software 
556    CALL DVRP_LOG_EVENT( 2, current_timestep_number )
557#endif
558
559    CALL location_message( 'starting timestep-sequence', .TRUE. )
560!
561!-- Start of the time loop
562    DO  WHILE ( simulated_time < end_time  .AND.  .NOT. stop_dt  .AND. &
563                .NOT. terminate_run )
564
565       CALL cpu_log( log_point_s(10), 'timesteps', 'start' )
566!
567!--    Vertical nesting: initialize fine grid
568       IF ( vnested ) THEN
569          IF ( .NOT. vnest_init  .AND.  simulated_time >= vnest_start_time )  THEN
570             CALL cpu_log( log_point(80), 'vnest_init', 'start' )
571             CALL vnest_init_fine
572             vnest_init = .TRUE.
573             CALL cpu_log( log_point(80), 'vnest_init', 'stop' )
574          ENDIF
575       ENDIF
576!
577!--    Determine ug, vg and w_subs in dependence on data from external file
578!--    LSF_DATA
579       IF ( large_scale_forcing .AND. lsf_vert )  THEN
580           CALL ls_forcing_vert ( simulated_time )
581           sums_ls_l = 0.0_wp
582       ENDIF
583
584!
585!--    Set pt_init and q_init to the current profiles taken from
586!--    NUDGING_DATA
587       IF ( nudging )  THEN
588           CALL nudge_ref ( simulated_time )
589!
590!--        Store temperature gradient at the top boundary for possible Neumann
591!--        boundary condition
592           bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1)
593           bc_q_t_val  = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1)
594       ENDIF
595!
596!--    If forcing by larger-scale models is applied, check if new data
597!--    at domain boundaries need to be read.
598       IF ( nesting_offline )  THEN
599          IF ( nest_offl%time(nest_offl%tind_p) <= time_since_reference_point )&
600             CALL netcdf_data_input_lsf
601       ENDIF
602
603!
604!--    Execute the gust module actions
605       IF ( gust_module_enabled )  THEN
606          CALL gust_actions( 'before_timestep' )
607       ENDIF
608
609!
610!--    Execute the user-defined actions
611       CALL user_actions( 'before_timestep' )
612
613!
614!--    Calculate forces by wind turbines
615       IF ( wind_turbine )  THEN
616
617          CALL cpu_log( log_point(55), 'wind_turbine', 'start' )
618
619          CALL wtm_forces
620
621          CALL cpu_log( log_point(55), 'wind_turbine', 'stop' )
622
623       ENDIF       
624       
625!
626!--    Start of intermediate step loop
627       intermediate_timestep_count = 0
628       DO  WHILE ( intermediate_timestep_count < &
629                   intermediate_timestep_count_max )
630
631          intermediate_timestep_count = intermediate_timestep_count + 1
632
633!
634!--       Set the steering factors for the prognostic equations which depend
635!--       on the timestep scheme
636          CALL timestep_scheme_steering
637
638!
639!--       Calculate those variables needed in the tendency terms which need
640!--       global communication
641          IF ( .NOT. use_single_reference_value  .AND. &
642               .NOT. use_initial_profile_as_reference )  THEN
643!
644!--          Horizontally averaged profiles to be used as reference state in
645!--          buoyancy terms (WARNING: only the respective last call of
646!--          calc_mean_profile defines the reference state!)
647             IF ( .NOT. neutral )  THEN
648                CALL calc_mean_profile( pt, 4 )
649                ref_state(:)  = hom(:,1,4,0) ! this is used in the buoyancy term
650             ENDIF
651             IF ( ocean )  THEN
652                CALL calc_mean_profile( rho_ocean, 64 )
653                ref_state(:)  = hom(:,1,64,0)
654             ENDIF
655             IF ( humidity )  THEN
656                CALL calc_mean_profile( vpt, 44 )
657                ref_state(:)  = hom(:,1,44,0)
658             ENDIF
659!
660!--          Assure that ref_state does not become zero at any level
661!--          ( might be the case if a vertical level is completely occupied
662!--            with topography ).
663             ref_state = MERGE( MAXVAL(ref_state), ref_state,                  &
664                                ref_state == 0.0_wp )
665          ENDIF
666
667          IF ( .NOT. constant_diffusion )  CALL production_e_init
668          IF ( ( ws_scheme_mom .OR. ws_scheme_sca )  .AND.  &
669               intermediate_timestep_count == 1 )  CALL ws_statistics
670!
671!--       In case of nudging calculate current nudging time scale and horizontal
672!--       means of u, v, pt and q
673          IF ( nudging )  THEN
674             CALL calc_tnudge( simulated_time )
675             CALL calc_mean_profile( u, 1 )
676             CALL calc_mean_profile( v, 2 )
677             CALL calc_mean_profile( pt, 4 )
678             CALL calc_mean_profile( q, 41 )
679          ENDIF
680
681!
682!--       Solve the prognostic equations. A fast cache optimized version with
683!--       only one single loop is used in case of Piascek-Williams advection
684!--       scheme. NEC vector machines use a different version, because
685!--       in the other versions a good vectorization is prohibited due to
686!--       inlining problems.
687          IF ( loop_optimization == 'cache' )  THEN
688             CALL prognostic_equations_cache
689          ELSEIF ( loop_optimization == 'vector' )  THEN
690             CALL prognostic_equations_vector
691          ENDIF
692
693!
694!--       Particle transport/physics with the Lagrangian particle model
695!--       (only once during intermediate steps, because it uses an Euler-step)
696!--       ### particle model should be moved before prognostic_equations, in order
697!--       to regard droplet interactions directly
698          IF ( particle_advection  .AND.                         &
699               simulated_time >= particle_advection_start  .AND. &
700               intermediate_timestep_count == 1 )  THEN
701             CALL lpm
702             first_call_lpm = .FALSE.
703          ENDIF
704
705!
706!--       Interaction of droplets with temperature and mixing ratio.
707!--       Droplet condensation and evaporation is calculated within
708!--       advec_particles.
709          IF ( cloud_droplets  .AND.  &
710               intermediate_timestep_count == intermediate_timestep_count_max )&
711          THEN
712             CALL interaction_droplets_ptq
713          ENDIF
714
715!
716!--       Movement of agents in multi agent system
717          IF ( agents_active  .AND.                                            &
718               simulated_time >= multi_agent_system_start  .AND.               &
719               intermediate_timestep_count == 1 )  THEN
720             CALL multi_agent_system
721             first_call_mas = .FALSE.
722          ENDIF
723
724!
725!--       Exchange of ghost points (lateral boundary conditions)
726          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'start' )
727
728          CALL exchange_horiz( u_p, nbgp )
729          CALL exchange_horiz( v_p, nbgp )
730          CALL exchange_horiz( w_p, nbgp )
731          CALL exchange_horiz( pt_p, nbgp )
732          IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e_p, nbgp )
733          IF ( rans_tke_e  .OR.  wang_kernel  .OR.  collision_turbulence       &
734               .OR.  use_sgs_for_particles )  THEN
735             IF ( rans_tke_e )  THEN
736                CALL exchange_horiz( diss_p, nbgp )
737             ELSE
738                CALL exchange_horiz( diss, nbgp )
739             ENDIF
740          ENDIF
741          IF ( ocean )  THEN
742             CALL exchange_horiz( sa_p, nbgp )
743             CALL exchange_horiz( rho_ocean, nbgp )
744             CALL exchange_horiz( prho, nbgp )
745          ENDIF
746          IF ( humidity )  THEN
747             CALL exchange_horiz( q_p, nbgp )
748             IF ( cloud_physics .AND. microphysics_morrison )  THEN
749                CALL exchange_horiz( qc_p, nbgp )
750                CALL exchange_horiz( nc_p, nbgp )
751             ENDIF
752             IF ( cloud_physics .AND. microphysics_seifert )  THEN
753                CALL exchange_horiz( qr_p, nbgp )
754                CALL exchange_horiz( nr_p, nbgp )
755             ENDIF
756          ENDIF
757          IF ( cloud_droplets )  THEN
758             CALL exchange_horiz( ql, nbgp )
759             CALL exchange_horiz( ql_c, nbgp )
760             CALL exchange_horiz( ql_v, nbgp )
761             CALL exchange_horiz( ql_vp, nbgp )
762          ENDIF
763          IF ( passive_scalar )  CALL exchange_horiz( s_p, nbgp )
764          IF ( air_chemistry )  THEN
765             DO  n = 1, nspec     
766                CALL exchange_horiz( chem_species(n)%conc_p, nbgp ) 
767             ENDDO
768          ENDIF
769
770          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' )
771
772!
773!--       Boundary conditions for the prognostic quantities (except of the
774!--       velocities at the outflow in case of a non-cyclic lateral wall)
775          CALL boundary_conds
776!
777!--       Swap the time levels in preparation for the next time step.
778          CALL swap_timelevel
779
780!
781!--       Vertical nesting: Interpolate fine grid data to the coarse grid
782          IF ( vnest_init ) THEN
783             CALL cpu_log( log_point(81), 'vnest_anterpolate', 'start' )
784             CALL vnest_anterpolate
785             CALL cpu_log( log_point(81), 'vnest_anterpolate', 'stop' )
786          ENDIF
787
788          IF ( nested_run )  THEN
789
790             CALL cpu_log( log_point(60), 'nesting', 'start' )
791!
792!--          Domain nesting. The data transfer subroutines pmci_parent_datatrans
793!--          and pmci_child_datatrans are called inside the wrapper
794!--          subroutine pmci_datatrans according to the control parameters
795!--          nesting_mode and nesting_datatransfer_mode.
796!--          TO_DO: why is nesting_mode given as a parameter here?
797             CALL pmci_datatrans( nesting_mode )
798
799             IF ( TRIM( nesting_mode ) == 'two-way' .OR.                       &
800                  nesting_mode == 'vertical' )  THEN
801!
802!--             Exchange_horiz is needed for all parent-domains after the
803!--             anterpolation
804                CALL exchange_horiz( u, nbgp )
805                CALL exchange_horiz( v, nbgp )
806                CALL exchange_horiz( w, nbgp )
807                IF ( .NOT. neutral )  CALL exchange_horiz( pt, nbgp )
808
809                IF ( humidity )  THEN
810
811                   CALL exchange_horiz( q, nbgp )
812
813                   IF ( cloud_physics  .AND.  microphysics_morrison )  THEN
814                       CALL exchange_horiz( qc, nbgp )
815                       CALL exchange_horiz( nc, nbgp )
816                   ENDIF
817                   IF ( cloud_physics  .AND.  microphysics_seifert )  THEN
818                       CALL exchange_horiz( qr, nbgp )
819                       CALL exchange_horiz( nr, nbgp )
820                   ENDIF
821
822                ENDIF
823
824                IF ( passive_scalar )  CALL exchange_horiz( s, nbgp )
825                IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e, nbgp )
826
827                IF ( .NOT. constant_diffusion  .AND.  rans_mode  .AND.         &
828                                                      rans_tke_e )             &
829                   CALL exchange_horiz( diss, nbgp )
830
831                IF ( air_chemistry )  THEN
832                   DO  n = 1, nspec     
833                      CALL exchange_horiz( chem_species(n)%conc, nbgp ) 
834                   ENDDO
835                ENDIF
836
837             ENDIF
838!
839!--          Set boundary conditions again after interpolation and anterpolation.
840             CALL pmci_boundary_conds
841!
842!--          Correct the w top-BC in nest domains to ensure mass conservation.
843!--          This action must never be done for the root domain. Vertical
844!--          nesting implies mass conservation.
845!--          Commented out April 18, 2018 as seemingly unnecessary.
846!--          Will later be completely removed.
847!--             IF ( child_domain )  THEN
848!--                CALL pmci_ensure_nest_mass_conservation
849!--             ENDIF
850
851
852             CALL cpu_log( log_point(60), 'nesting', 'stop' )
853
854          ENDIF
855
856!
857!--       Temperature offset must be imposed at cyclic boundaries in x-direction
858!--       when a sloping surface is used
859          IF ( sloping_surface )  THEN
860             IF ( nxl ==  0 )  pt(:,:,nxlg:nxl-1) = pt(:,:,nxlg:nxl-1) - &
861                                                    pt_slope_offset
862             IF ( nxr == nx )  pt(:,:,nxr+1:nxrg) = pt(:,:,nxr+1:nxrg) + &
863                                                    pt_slope_offset
864          ENDIF
865
866!
867!--       Impose a turbulent inflow using the recycling method
868          IF ( turbulent_inflow )  CALL  inflow_turbulence
869
870!
871!--       Impose a turbulent inflow using synthetic generated turbulence
872          IF ( use_syn_turb_gen )  CALL  stg_main
873
874!
875!--       Set values at outflow boundary using the special outflow condition
876          IF ( turbulent_outflow )  CALL  outflow_turbulence
877
878!
879!--       Impose a random perturbation on the horizontal velocity field
880          IF ( create_disturbances  .AND.                                      &
881               ( call_psolver_at_all_substeps  .AND.                           &
882               intermediate_timestep_count == intermediate_timestep_count_max )&
883          .OR. ( .NOT. call_psolver_at_all_substeps  .AND.                     &
884               intermediate_timestep_count == 1 ) )                            &
885          THEN
886             time_disturb = time_disturb + dt_3d
887             IF ( time_disturb >= dt_disturb )  THEN
888                IF ( disturbance_energy_limit /= 0.0_wp  .AND.                 &
889                     hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit )  THEN
890                   CALL disturb_field( 'u', tend, u )
891                   CALL disturb_field( 'v', tend, v )
892                ELSEIF ( ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )            &
893                     .AND. .NOT. child_domain  .AND.  .NOT.  nesting_offline )  &
894                THEN
895!
896!--                Runs with a non-cyclic lateral wall need perturbations
897!--                near the inflow throughout the whole simulation
898                   dist_range = 1
899                   CALL disturb_field( 'u', tend, u )
900                   CALL disturb_field( 'v', tend, v )
901                   dist_range = 0
902                ENDIF
903                time_disturb = time_disturb - dt_disturb
904             ENDIF
905          ENDIF
906
907!
908!--       Map forcing data derived from larger scale model onto domain
909!--       boundaries.
910          IF ( nesting_offline  .AND.  intermediate_timestep_count ==          &
911                                       intermediate_timestep_count_max  )  THEN
912             CALL lsf_nesting_offline
913!
914!--          Moreover, ensure mass conservation
915             CALL lsf_nesting_offline_mass_conservation
916          ENDIF
917
918!
919!--       Reduce the velocity divergence via the equation for perturbation
920!--       pressure.
921          IF ( intermediate_timestep_count == 1  .OR. &
922                call_psolver_at_all_substeps )  THEN
923
924             IF (  vnest_init ) THEN
925!
926!--             Compute pressure in the CG, interpolate top boundary conditions
927!--             to the FG and then compute pressure in the FG
928                IF ( coupling_mode == 'vnested_crse' )  CALL pres
929
930                CALL cpu_log( log_point(82), 'vnest_bc', 'start' )
931                CALL vnest_boundary_conds
932                CALL cpu_log( log_point(82), 'vnest_bc', 'stop' )
933 
934                IF ( coupling_mode == 'vnested_fine' )  CALL pres
935
936!--             Anterpolate TKE, satisfy Germano Identity
937                CALL cpu_log( log_point(83), 'vnest_anter_e', 'start' )
938                CALL vnest_anterpolate_e
939                CALL cpu_log( log_point(83), 'vnest_anter_e', 'stop' )
940
941             ELSE
942
943                CALL pres
944
945             ENDIF
946
947          ENDIF
948
949!
950!--       If required, compute liquid water content
951          IF ( cloud_physics )  THEN
952             CALL calc_liquid_water_content
953          ENDIF
954!
955!--       If required, compute virtual potential temperature
956          IF ( humidity )  THEN
957             CALL compute_vpt 
958          ENDIF 
959
960!
961!--       Compute the diffusion quantities
962          IF ( .NOT. constant_diffusion )  THEN
963
964!
965!--          Determine surface fluxes shf and qsws and surface values
966!--          pt_surface and q_surface in dependence on data from external
967!--          file LSF_DATA respectively
968             IF ( ( large_scale_forcing .AND. lsf_surf ) .AND. &
969                 intermediate_timestep_count == intermediate_timestep_count_max )&
970             THEN
971                CALL ls_forcing_surf( simulated_time )
972             ENDIF
973
974!
975!--          First the vertical (and horizontal) fluxes in the surface
976!--          (constant flux) layer are computed
977             IF ( constant_flux_layer )  THEN
978                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'start' )
979                CALL surface_layer_fluxes
980                CALL cpu_log( log_point(19), 'surface_layer_fluxes', 'stop' )
981             ENDIF
982!
983!--          If required, solve the energy balance for the surface and run soil
984!--          model. Call for horizontal as well as vertical surfaces
985             IF ( land_surface .AND. time_since_reference_point >= skip_time_do_lsm)  THEN
986
987                CALL cpu_log( log_point(54), 'land_surface', 'start' )
988!
989!--             Call for horizontal upward-facing surfaces
990                CALL lsm_energy_balance( .TRUE., -1 )
991                CALL lsm_soil_model( .TRUE., -1, .TRUE. )
992!
993!--             Call for northward-facing surfaces
994                CALL lsm_energy_balance( .FALSE., 0 )
995                CALL lsm_soil_model( .FALSE., 0, .TRUE. )
996!
997!--             Call for southward-facing surfaces
998                CALL lsm_energy_balance( .FALSE., 1 )
999                CALL lsm_soil_model( .FALSE., 1, .TRUE. )
1000!
1001!--             Call for eastward-facing surfaces
1002                CALL lsm_energy_balance( .FALSE., 2 )
1003                CALL lsm_soil_model( .FALSE., 2, .TRUE. )
1004!
1005!--             Call for westward-facing surfaces
1006                CALL lsm_energy_balance( .FALSE., 3 )
1007                CALL lsm_soil_model( .FALSE., 3, .TRUE. )
1008!
1009!--             At the end, set boundary conditons for potential temperature
1010!--             and humidity after running the land-surface model. This
1011!--             might be important for the nesting, where arrays are transfered.
1012                CALL lsm_boundary_condition
1013
1014                CALL cpu_log( log_point(54), 'land_surface', 'stop' )
1015             ENDIF
1016!
1017!--          If required, solve the energy balance for urban surfaces and run
1018!--          the material heat model
1019             IF (urban_surface) THEN
1020                CALL cpu_log( log_point(74), 'urban_surface', 'start' )
1021               
1022                CALL usm_surface_energy_balance
1023                IF ( usm_material_model )  THEN
1024                   CALL usm_green_heat_model
1025                   CALL usm_material_heat_model
1026                ENDIF
1027
1028                CALL usm_temperature_near_surface
1029!
1030!--             At the end, set boundary conditons for potential temperature
1031!--             and humidity after running the urban-surface model. This
1032!--             might be important for the nesting, where arrays are transfered.
1033                CALL usm_boundary_condition
1034
1035                CALL cpu_log( log_point(74), 'urban_surface', 'stop' )
1036             ENDIF
1037!
1038!--          Compute the diffusion coefficients
1039             CALL cpu_log( log_point(17), 'diffusivities', 'start' )
1040             IF ( .NOT. humidity ) THEN
1041                IF ( ocean )  THEN
1042                   CALL tcm_diffusivities( prho, prho_reference )
1043                ELSE
1044                   CALL tcm_diffusivities( pt, pt_reference )
1045                ENDIF
1046             ELSE
1047                CALL tcm_diffusivities( vpt, pt_reference )
1048             ENDIF
1049             CALL cpu_log( log_point(17), 'diffusivities', 'stop' )
1050!
1051!--          Vertical nesting: set fine grid eddy viscosity top boundary condition
1052             IF ( vnest_init )  CALL vnest_boundary_conds_khkm
1053
1054          ENDIF
1055
1056!
1057!--       If required, calculate radiative fluxes and heating rates
1058          IF ( radiation .AND. intermediate_timestep_count                     &
1059               == intermediate_timestep_count_max .AND. time_since_reference_point >    &
1060               skip_time_do_radiation )  THEN
1061
1062               time_radiation = time_radiation + dt_3d
1063
1064             IF ( time_radiation >= dt_radiation .OR. force_radiation_call )   &
1065             THEN
1066
1067                CALL cpu_log( log_point(50), 'radiation', 'start' )
1068
1069                IF ( .NOT. force_radiation_call )  THEN
1070                   time_radiation = time_radiation - dt_radiation
1071                ENDIF
1072
1073!
1074!--             Adjust the current_ time to the time step of the radiation model.
1075!--             Needed since radiation is pre-calculated and stored only on apparent
1076!--             solar positions
1077                it = FLOOR(time_since_reference_point/dt_radiation)
1078                tsrp_org = time_since_reference_point
1079                time_since_reference_point = REAL(it,wp) * dt_radiation
1080
1081                CALL radiation_control
1082
1083                CALL cpu_log( log_point(50), 'radiation', 'stop' )
1084
1085                IF ( urban_surface  .OR.  land_surface  .AND.                  &
1086                     radiation_interactions )  THEN
1087                   CALL cpu_log( log_point(75), 'radiation_interaction', 'start' )
1088                   CALL radiation_interaction
1089                   CALL cpu_log( log_point(75), 'radiation_interaction', 'stop' )
1090                ENDIF
1091   
1092!
1093!--             Return the current time to its original value
1094                time_since_reference_point = tsrp_org
1095
1096             ENDIF
1097          ENDIF
1098
1099       ENDDO   ! Intermediate step loop
1100!
1101!--    If required, consider chemical emissions
1102!--    (todo (FK): Implement hourly call of emissions, using time_utc from
1103!--                data_and_time_mod.f90;
1104!--                move the CALL to appropriate location)
1105       IF ( air_chemistry ) THEN
1106          CALL chem_emissions
1107       ENDIF
1108!
1109!--    If required, do UV exposure calculations
1110       IF ( uv_exposure )  THEN
1111          CALL uvem_calc_exposure
1112       ENDIF
1113!
1114!--    Increase simulation time and output times
1115       nr_timesteps_this_run      = nr_timesteps_this_run + 1
1116       current_timestep_number    = current_timestep_number + 1
1117       simulated_time             = simulated_time   + dt_3d
1118       time_since_reference_point = simulated_time - coupling_start_time
1119       simulated_time_chr         = time_to_string( time_since_reference_point )
1120
1121
1122
1123
1124       IF ( simulated_time >= skip_time_data_output_av )  THEN
1125          time_do_av         = time_do_av       + dt_3d
1126       ENDIF
1127       IF ( simulated_time >= skip_time_do2d_xy )  THEN
1128          time_do2d_xy       = time_do2d_xy     + dt_3d
1129       ENDIF
1130       IF ( simulated_time >= skip_time_do2d_xz )  THEN
1131          time_do2d_xz       = time_do2d_xz     + dt_3d
1132       ENDIF
1133       IF ( simulated_time >= skip_time_do2d_yz )  THEN
1134          time_do2d_yz       = time_do2d_yz     + dt_3d
1135       ENDIF
1136       IF ( simulated_time >= skip_time_do3d    )  THEN
1137          time_do3d          = time_do3d        + dt_3d
1138       ENDIF
1139       DO  mid = 1, masks
1140          IF ( simulated_time >= skip_time_domask(mid) )  THEN
1141             time_domask(mid)= time_domask(mid) + dt_3d
1142          ENDIF
1143       ENDDO
1144       time_dvrp          = time_dvrp        + dt_3d
1145       IF ( simulated_time >= skip_time_dosp )  THEN
1146          time_dosp       = time_dosp        + dt_3d
1147       ENDIF
1148       time_dots          = time_dots        + dt_3d
1149       IF ( .NOT. first_call_lpm )  THEN
1150          time_dopts      = time_dopts       + dt_3d
1151       ENDIF
1152       IF ( simulated_time >= skip_time_dopr )  THEN
1153          time_dopr       = time_dopr        + dt_3d
1154       ENDIF
1155       time_dopr_listing          = time_dopr_listing        + dt_3d
1156       time_run_control   = time_run_control + dt_3d
1157
1158!
1159!--    Data exchange between coupled models
1160       IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled                   &
1161                                          .AND. .NOT. vnested )  THEN
1162          time_coupling = time_coupling + dt_3d
1163
1164!
1165!--       In case of model termination initiated by the local model
1166!--       (terminate_coupled > 0), the coupler must be skipped because it would
1167!--       cause an MPI intercomminucation hang.
1168!--       If necessary, the coupler will be called at the beginning of the
1169!--       next restart run.
1170          DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 )
1171             CALL surface_coupler
1172             time_coupling = time_coupling - dt_coupling
1173          ENDDO
1174       ENDIF
1175
1176!
1177!--    Execute the gust module actions
1178       IF ( gust_module_enabled )  THEN
1179          CALL gust_actions( 'after_integration' )
1180       ENDIF
1181
1182!
1183!--    Execute user-defined actions
1184       CALL user_actions( 'after_integration' )
1185
1186!
1187!--    If Galilei transformation is used, determine the distance that the
1188!--    model has moved so far
1189       IF ( galilei_transformation )  THEN
1190          advected_distance_x = advected_distance_x + u_gtrans * dt_3d
1191          advected_distance_y = advected_distance_y + v_gtrans * dt_3d
1192       ENDIF
1193
1194!
1195!--    Check, if restart is necessary (because cpu-time is expiring or
1196!--    because it is forced by user) and set stop flag
1197!--    This call is skipped if the remote model has already initiated a restart.
1198       IF ( .NOT. terminate_run )  CALL check_for_restart
1199
1200!
1201!--    Carry out statistical analysis and output at the requested output times.
1202!--    The MOD function is used for calculating the output time counters (like
1203!--    time_dopr) in order to regard a possible decrease of the output time
1204!--    interval in case of restart runs
1205
1206!
1207!--    Set a flag indicating that so far no statistics have been created
1208!--    for this time step
1209       flow_statistics_called = .FALSE.
1210
1211!
1212!--    If required, call flow_statistics for averaging in time
1213       IF ( averaging_interval_pr /= 0.0_wp  .AND.  &
1214            ( dt_dopr - time_dopr ) <= averaging_interval_pr  .AND.  &
1215            simulated_time >= skip_time_dopr )  THEN
1216          time_dopr_av = time_dopr_av + dt_3d
1217          IF ( time_dopr_av >= dt_averaging_input_pr )  THEN
1218             do_sum = .TRUE.
1219             time_dopr_av = MOD( time_dopr_av, &
1220                                    MAX( dt_averaging_input_pr, dt_3d ) )
1221          ENDIF
1222       ENDIF
1223       IF ( do_sum )  CALL flow_statistics
1224
1225!
1226!--    Sum-up 3d-arrays for later output of time-averaged 2d/3d/masked data
1227       IF ( averaging_interval /= 0.0_wp  .AND.                                &
1228            ( dt_data_output_av - time_do_av ) <= averaging_interval  .AND. &
1229            simulated_time >= skip_time_data_output_av )                    &
1230       THEN
1231          time_do_sla = time_do_sla + dt_3d
1232          IF ( time_do_sla >= dt_averaging_input )  THEN
1233             CALL sum_up_3d_data
1234             average_count_3d = average_count_3d + 1
1235             time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) )
1236          ENDIF
1237       ENDIF
1238
1239!
1240!--    Calculate spectra for time averaging
1241       IF ( averaging_interval_sp /= 0.0_wp  .AND.  &
1242            ( dt_dosp - time_dosp ) <= averaging_interval_sp  .AND.  &
1243            simulated_time >= skip_time_dosp )  THEN
1244          time_dosp_av = time_dosp_av + dt_3d
1245          IF ( time_dosp_av >= dt_averaging_input_pr )  THEN
1246             CALL calc_spectra
1247             time_dosp_av = MOD( time_dosp_av, &
1248                                 MAX( dt_averaging_input_pr, dt_3d ) )
1249          ENDIF
1250       ENDIF
1251
1252!
1253!--    Call flight module and output data
1254       IF ( virtual_flight )  THEN
1255          CALL flight_measurement
1256          CALL data_output_flight
1257       ENDIF
1258
1259!
1260!--    Profile output (ASCII) on file
1261       IF ( time_dopr_listing >= dt_dopr_listing )  THEN
1262          CALL print_1d
1263          time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, &
1264                                                           dt_3d ) )
1265       ENDIF
1266
1267!
1268!--    Graphic output for PROFIL
1269       IF ( time_dopr >= dt_dopr )  THEN
1270          IF ( dopr_n /= 0 )  CALL data_output_profiles
1271          time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) )
1272          time_dopr_av = 0.0_wp    ! due to averaging (see above)
1273       ENDIF
1274
1275!
1276!--    Graphic output for time series
1277       IF ( time_dots >= dt_dots )  THEN
1278          CALL data_output_tseries
1279          time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) )
1280       ENDIF
1281
1282!
1283!--    Output of spectra (formatted for use with PROFIL), in case of no
1284!--    time averaging, spectra has to be calculated before
1285       IF ( time_dosp >= dt_dosp )  THEN
1286          IF ( average_count_sp == 0 )  CALL calc_spectra
1287          CALL data_output_spectra
1288          time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) )
1289       ENDIF
1290
1291!
1292!--    2d-data output (cross-sections)
1293       IF ( time_do2d_xy >= dt_do2d_xy )  THEN
1294          CALL data_output_2d( 'xy', 0 )
1295          time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) )
1296       ENDIF
1297       IF ( time_do2d_xz >= dt_do2d_xz )  THEN
1298          CALL data_output_2d( 'xz', 0 )
1299          time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) )
1300       ENDIF
1301       IF ( time_do2d_yz >= dt_do2d_yz )  THEN
1302          CALL data_output_2d( 'yz', 0 )
1303          time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) )
1304       ENDIF
1305
1306!
1307!--    3d-data output (volume data)
1308       IF ( time_do3d >= dt_do3d )  THEN
1309          CALL data_output_3d( 0 )
1310          time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) )
1311       ENDIF
1312
1313!
1314!--    Masked data output
1315       DO  mid = 1, masks
1316          IF ( time_domask(mid) >= dt_domask(mid) )  THEN
1317             CALL data_output_mask( 0 )
1318             time_domask(mid) = MOD( time_domask(mid),  &
1319                                     MAX( dt_domask(mid), dt_3d ) )
1320          ENDIF
1321       ENDDO
1322
1323!
1324!--    Output of time-averaged 2d/3d/masked data
1325       IF ( time_do_av >= dt_data_output_av )  THEN
1326          CALL average_3d_data
1327          CALL data_output_2d( 'xy', 1 )
1328          CALL data_output_2d( 'xz', 1 )
1329          CALL data_output_2d( 'yz', 1 )
1330          CALL data_output_3d( 1 )
1331          DO  mid = 1, masks
1332             CALL data_output_mask( 1 )
1333          ENDDO
1334          time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) )
1335       ENDIF
1336
1337!
1338!--    Output of particle time series
1339       IF ( particle_advection )  THEN
1340          IF ( time_dopts >= dt_dopts  .OR. &
1341               ( simulated_time >= particle_advection_start  .AND. &
1342                 first_call_lpm ) )  THEN
1343             CALL data_output_ptseries
1344             time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
1345          ENDIF
1346       ENDIF
1347
1348!
1349!--    Output of dvrp-graphics (isosurface, particles, slicer)
1350#if defined( __dvrp_graphics )
1351       CALL DVRP_LOG_EVENT( -2, current_timestep_number-1 )
1352#endif
1353       IF ( time_dvrp >= dt_dvrp )  THEN
1354          CALL data_output_dvrp
1355          time_dvrp = MOD( time_dvrp, MAX( dt_dvrp, dt_3d ) )
1356       ENDIF
1357#if defined( __dvrp_graphics )
1358       CALL DVRP_LOG_EVENT( 2, current_timestep_number )
1359#endif
1360
1361!
1362!--    If required, set the heat flux for the next time step at a random value
1363       IF ( constant_heatflux  .AND.  random_heatflux )  THEN
1364          IF ( surf_def_h(0)%ns >= 1 )  CALL disturb_heatflux( surf_def_h(0) )
1365          IF ( surf_lsm_h%ns    >= 1 )  CALL disturb_heatflux( surf_lsm_h    )
1366          IF ( surf_usm_h%ns    >= 1 )  CALL disturb_heatflux( surf_usm_h    )
1367       ENDIF
1368
1369!
1370!--    Execute the gust module actions
1371       IF ( gust_module_enabled )  THEN
1372          CALL gust_actions( 'after_timestep' )
1373       ENDIF
1374
1375!
1376!--    Execute user-defined actions
1377       CALL user_actions( 'after_timestep' )
1378
1379!
1380!--    Determine size of next time step. Save timestep dt_3d because it is
1381!--    newly calculated in routine timestep, but required further below for
1382!--    steering the run control output interval
1383       dt_3d_old = dt_3d
1384       CALL timestep
1385
1386!
1387!--    Synchronize the timestep in case of nested run.
1388       IF ( nested_run )  THEN
1389!
1390!--       Synchronize by unifying the time step.
1391!--       Global minimum of all time-steps is used for all.
1392          CALL pmci_synchronize
1393       ENDIF
1394
1395!
1396!--    Computation and output of run control parameters.
1397!--    This is also done whenever perturbations have been imposed
1398       IF ( time_run_control >= dt_run_control  .OR.                     &
1399            timestep_scheme(1:5) /= 'runge'  .OR.  disturbance_created ) &
1400       THEN
1401          CALL run_control
1402          IF ( time_run_control >= dt_run_control )  THEN
1403             time_run_control = MOD( time_run_control, &
1404                                     MAX( dt_run_control, dt_3d_old ) )
1405          ENDIF
1406       ENDIF
1407
1408!
1409!--    Output elapsed simulated time in form of a progress bar on stdout
1410       IF ( myid == 0 )  CALL output_progress_bar
1411
1412       CALL cpu_log( log_point_s(10), 'timesteps', 'stop' )
1413
1414
1415    ENDDO   ! time loop
1416
1417!-- Vertical nesting: Deallocate variables initialized for vertical nesting   
1418    IF ( vnest_init )  CALL vnest_deallocate
1419
1420    IF ( myid == 0 )  CALL finish_progress_bar
1421
1422#if defined( __dvrp_graphics )
1423    CALL DVRP_LOG_EVENT( -2, current_timestep_number )
1424#endif
1425
1426    CALL location_message( 'finished time-stepping', .TRUE. )
1427
1428 END SUBROUTINE time_integration
Note: See TracBrowser for help on using the repository browser.