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

Last change on this file since 2801 was 2801, checked in by thiele, 6 years ago

Introduce particle transfer in nested models

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