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

Last change on this file since 2303 was 2299, checked in by maronga, 7 years ago

improvements for spinup mechanism

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