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

Last change on this file since 2543 was 2365, checked in by kanani, 7 years ago

Vertical nesting implemented (SadiqHuq?)

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