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

Last change on this file since 2605 was 2563, checked in by Giersch, 7 years ago

Restart runs with the usage of the wind turbine model are possible now. Further small at reading/writing restart data

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