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

Last change on this file since 2769 was 2766, checked in by kanani, 7 years ago

Removal of chem directive, plus minor changes

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