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

Last change on this file since 2636 was 2617, checked in by suehring, 7 years ago

Avoid that the mean reference state in buoyancy term does not become zero

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