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

Last change on this file since 1919 was 1919, checked in by raasch, 8 years ago

last commit documented

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