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

Last change on this file since 2742 was 2718, checked in by maronga, 7 years ago

deleting of deprecated files; headers updated where needed

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