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

Last change on this file since 3039 was 3014, checked in by maronga, 7 years ago

series of bugfixes

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