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

Last change on this file since 1036 was 1036, checked in by raasch, 12 years ago

code has been put under the GNU General Public License (v3)

  • Property svn:keywords set to Id
File size: 21.9 KB
RevLine 
[1]1 SUBROUTINE time_integration
2
[1036]3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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-2012  Leibniz University Hannover
18!--------------------------------------------------------------------------------!
19!
[484]20! Current revisions:
[1]21! -----------------
[392]22!
[1020]23!
[392]24! Former revisions:
25! -----------------
26! $Id: time_integration.f90 1036 2012-10-22 13:43:42Z raasch $
27!
[1020]28! 1019 2012-09-28 06:46:45Z raasch
29! non-optimized version of prognostic_equations removed
30!
[1017]31! 1015 2012-09-27 09:23:24Z raasch
32! +call of prognostic_equations_acc
33!
[1002]34! 1001 2012-09-13 14:08:46Z raasch
35! all actions concerning leapfrog- and upstream-spline-scheme removed
36!
[850]37! 849 2012-03-15 10:35:09Z raasch
38! advec_particles renamed lpm, first_call_advec_particles renamed first_call_lpm
39!
[826]40! 825 2012-02-19 03:03:44Z raasch
41! wang_collision_kernel renamed wang_kernel
42!
[791]43! 790 2011-11-29 03:11:20Z raasch
44! exchange of ghostpoints for array diss
45!
[708]46! 707 2011-03-29 11:39:40Z raasch
47! bc_lr/ns replaced by bc_lr/ns_cyc, calls of exchange_horiz are modified,
48! adaption to sloping surface
49!
[668]50! 667  2010-12-23 12:06:00Z suehring/gryschka
51! Calls of exchange_horiz are modified.
52! Adaption to slooping surface.
53!
[482]54! 449 2010-02-02 11:23:59Z raasch
55! Bugfix: exchange of ghost points for prho included
56!
[449]57! 410 2009-12-04 17:05:40Z letzel
58! masked data output
59!
[392]60! 388 2009-09-23 09:40:33Z raasch
[388]61! Using prho instead of rho in diffusvities.
[291]62! Coupling with independent precursor runs.
63! Bugfix: output of particle time series only if particle advection is switched
[253]64!         on
[110]65!
[198]66! 151 2008-03-07 13:42:18Z raasch
67! inflow turbulence is imposed by calling new routine inflow_turbulence
68!
[110]69! 108 2007-08-24 15:10:38Z letzel
[106]70! Call of new routine surface_coupler,
71! presure solver is called after the first Runge-Kutta substep instead of the
72! last in case that call_psolver_at_all_substeps = .F.; for this case, the
73! random perturbation has to be added to the velocity fields also after the
74! first substep
[77]75!
[98]76! 97 2007-06-21 08:23:15Z raasch
77! diffusivities is called with argument rho in case of ocean runs,
78! new argument pt_/prho_reference in calls of diffusivities,
79! ghostpoint exchange for salinity and density
80!
[90]81! 87 2007-05-22 15:46:47Z raasch
82! var_hom renamed pr_palm
83!
[77]84! 75 2007-03-22 09:54:05Z raasch
[46]85! Move call of user_actions( 'after_integration' ) below increment of times
[63]86! and counters,
87! calls of prognostic_equations_.. changed to .._noopt, .._cache, and
[75]88! .._vector, these calls are now controlled by switch loop_optimization,
89! uxrp, vynp eliminated, 2nd+3rd argument removed from exchange horiz,
90! moisture renamed humidity
[1]91!
[3]92! RCS Log replace by Id keyword, revision history cleaned up
93!
[1]94! Revision 1.8  2006/08/22 14:16:05  raasch
95! Disturbances are imposed only for the last Runge-Kutta-substep
96!
97! Revision 1.2  2004/04/30 13:03:40  raasch
98! decalpha-specific warning removed, routine name changed to time_integration,
99! particle advection is carried out only once during the intermediate steps,
100! impulse_advec renamed momentum_advec
101!
102! Revision 1.1  1997/08/11 06:19:04  raasch
103! Initial revision
104!
105!
106! Description:
107! ------------
108! Integration in time of the model equations, statistical analysis and graphic
109! output
110!------------------------------------------------------------------------------!
111
112    USE arrays_3d
113    USE averaging
114    USE control_parameters
115    USE cpulog
116#if defined( __dvrp_graphics )
117    USE DVRP
118#endif
119    USE grid_variables
120    USE indices
121    USE interaction_droplets_ptq_mod
122    USE interfaces
123    USE particle_attributes
124    USE pegrid
125    USE prognostic_equations_mod
126    USE statistics
127    USE user_actions_mod
128
129    IMPLICIT NONE
130
131    CHARACTER (LEN=9) ::  time_to_string
132    INTEGER ::  i, j, k
133
134!
135!-- At the beginning of a simulation determine the time step as well as
136!-- determine and print out the run control parameters
137    IF ( simulated_time == 0.0 )  CALL timestep
[667]138
[1]139    CALL run_control
140
[667]141
[108]142!
143!-- Data exchange between coupled models in case that a call has been omitted
144!-- at the end of the previous run of a job chain.
[291]145    IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
[108]146!
147!--    In case of model termination initiated by the local model the coupler
148!--    must not be called because this would again cause an MPI hang.
149       DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 )
150          CALL surface_coupler
151          time_coupling = time_coupling - dt_coupling
152       ENDDO
[348]153       IF (time_coupling == 0.0 .AND. time_since_reference_point < dt_coupling)&
154       THEN
155          time_coupling = time_since_reference_point
156       ENDIF
[108]157    ENDIF
158
159
[1]160#if defined( __dvrp_graphics )
161!
162!-- Time measurement with dvrp software 
163    CALL DVRP_LOG_EVENT( 2, current_timestep_number )
164#endif
165
166!
167!-- Start of the time loop
168    DO  WHILE ( simulated_time < end_time  .AND.  .NOT. stop_dt  .AND. &
169                .NOT. terminate_run )
170
171       CALL cpu_log( log_point_s(10), 'timesteps', 'start' )
172!
173!--    Determine size of next time step
174       IF ( simulated_time /= 0.0 )  CALL timestep
175!
176!--    Execute the user-defined actions
177       CALL user_actions( 'before_timestep' )
178
179!
180!--    Start of intermediate step loop
181       intermediate_timestep_count = 0
182       DO  WHILE ( intermediate_timestep_count < &
183                   intermediate_timestep_count_max )
184
185          intermediate_timestep_count = intermediate_timestep_count + 1
186
187!
188!--       Set the steering factors for the prognostic equations which depend
189!--       on the timestep scheme
190          CALL timestep_scheme_steering
191
192!
193!--       Solve the prognostic equations. A fast cache optimized version with
194!--       only one single loop is used in case of Piascek-Williams advection
195!--       scheme. NEC vector machines use a different version, because
196!--       in the other versions a good vectorization is prohibited due to
197!--       inlining problems.
[1019]198          IF ( loop_optimization == 'cache' )  THEN
199             CALL prognostic_equations_cache
200          ELSEIF ( loop_optimization == 'vector' )  THEN
[63]201             CALL prognostic_equations_vector
[1015]202          ELSEIF ( loop_optimization == 'acc' )  THEN
203             CALL prognostic_equations_acc
[1]204          ENDIF
205
206!
[849]207!--       Particle transport/physics with the Lagrangian particle model
208!--       (only once during intermediate steps, because it uses an Euler-step)
[63]209          IF ( particle_advection  .AND.                         &
210               simulated_time >= particle_advection_start  .AND. &
[1]211               intermediate_timestep_count == 1 )  THEN
[849]212             CALL lpm
213             first_call_lpm = .FALSE.
[1]214          ENDIF
215
216!
217!--       Interaction of droplets with temperature and specific humidity.
218!--       Droplet condensation and evaporation is calculated within
219!--       advec_particles.
220          IF ( cloud_droplets  .AND.  &
221               intermediate_timestep_count == intermediate_timestep_count_max )&
222          THEN
223             CALL interaction_droplets_ptq
224          ENDIF
225
226!
227!--       Exchange of ghost points (lateral boundary conditions)
228          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'start' )
[667]229          CALL exchange_horiz( u_p, nbgp )
230          CALL exchange_horiz( v_p, nbgp )
231          CALL exchange_horiz( w_p, nbgp )
232          CALL exchange_horiz( pt_p, nbgp )
233          IF ( .NOT. constant_diffusion )  CALL exchange_horiz( e_p, nbgp )
[95]234          IF ( ocean )  THEN
[667]235             CALL exchange_horiz( sa_p, nbgp )
236             CALL exchange_horiz( rho, nbgp )
237             CALL exchange_horiz( prho, nbgp )
[95]238          ENDIF
[790]239          IF (humidity  .OR.  passive_scalar)  CALL exchange_horiz( q_p, nbgp )
[1]240          IF ( cloud_droplets )  THEN
[667]241             CALL exchange_horiz( ql, nbgp )
242             CALL exchange_horiz( ql_c, nbgp )
243             CALL exchange_horiz( ql_v, nbgp )
244             CALL exchange_horiz( ql_vp, nbgp )
[1]245          ENDIF
[825]246          IF ( wang_kernel )  CALL exchange_horiz( diss, nbgp )
[1]247
248          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' )
249
250!
251!--       Boundary conditions for the prognostic quantities (except of the
252!--       velocities at the outflow in case of a non-cyclic lateral wall)
253          CALL boundary_conds( 'main' )
254
255!
[73]256!--       Swap the time levels in preparation for the next time step.
257          CALL swap_timelevel
258
259!
[1]260!--       Temperature offset must be imposed at cyclic boundaries in x-direction
261!--       when a sloping surface is used
262          IF ( sloping_surface )  THEN
[707]263             IF ( nxl ==  0 )  pt(:,:,nxlg:nxl-1) = pt(:,:,nxlg:nxl-1) - &
264                                                    pt_slope_offset
265             IF ( nxr == nx )  pt(:,:,nxr+1:nxrg) = pt(:,:,nxr+1:nxrg) + &
266                                                    pt_slope_offset
[1]267          ENDIF
268
269!
[151]270!--       Impose a turbulent inflow using the recycling method
271          IF ( turbulent_inflow )  CALL  inflow_turbulence
272
273!
[1]274!--       Impose a random perturbation on the horizontal velocity field
[106]275          IF ( create_disturbances  .AND.                                      &
276               ( call_psolver_at_all_substeps  .AND.                           &
[1]277               intermediate_timestep_count == intermediate_timestep_count_max )&
[106]278          .OR. ( .NOT. call_psolver_at_all_substeps  .AND.                     &
279               intermediate_timestep_count == 1 ) )                            &
[1]280          THEN
281             time_disturb = time_disturb + dt_3d
282             IF ( time_disturb >= dt_disturb )  THEN
[87]283                IF ( hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit )  THEN
[75]284                   CALL disturb_field( nzb_u_inner, tend, u )
285                   CALL disturb_field( nzb_v_inner, tend, v )
[707]286                ELSEIF ( .NOT. bc_lr_cyc  .OR.  .NOT. bc_ns_cyc )  THEN
[1]287!
288!--                Runs with a non-cyclic lateral wall need perturbations
289!--                near the inflow throughout the whole simulation
290                   dist_range = 1
[75]291                   CALL disturb_field( nzb_u_inner, tend, u )
292                   CALL disturb_field( nzb_v_inner, tend, v )
[1]293                   dist_range = 0
294                ENDIF
295                time_disturb = time_disturb - dt_disturb
296             ENDIF
297          ENDIF
298
299!
300!--       Reduce the velocity divergence via the equation for perturbation
301!--       pressure.
[106]302          IF ( intermediate_timestep_count == 1  .OR. &
303                call_psolver_at_all_substeps )  THEN
[1]304             CALL pres
305          ENDIF
[1015]306!
307!--       Update device memory for calculating diffusion quantities and for next
308!--       timestep
309          !$acc update device( e, pt, u, v, w )
310          !$acc update device( q )  if ( allocated( q ) )
[1]311
312!
313!--       If required, compute virtuell potential temperature
[1015]314          IF ( humidity )  THEN
315             CALL compute_vpt
316             !$acc update device( vpt )
317          ENDIF
[1]318
319!
320!--       If required, compute liquid water content
[1015]321          IF ( cloud_physics )  THEN
322             CALL calc_liquid_water_content
323             !$acc update device( ql )
324          ENDIF
[1]325
326!
327!--       Compute the diffusion quantities
328          IF ( .NOT. constant_diffusion )  THEN
329
330!
331!--          First the vertical fluxes in the Prandtl layer are being computed
332             IF ( prandtl_layer )  THEN
333                CALL cpu_log( log_point(19), 'prandtl_fluxes', 'start' )
334                CALL prandtl_fluxes
335                CALL cpu_log( log_point(19), 'prandtl_fluxes', 'stop' )
[1015]336!
337!++             Statistics still require updates on host
338                !$acc update host( qs, qsws, rif, shf, ts )
[1]339             ENDIF
340
341!
342!--          Compute the diffusion coefficients
343             CALL cpu_log( log_point(17), 'diffusivities', 'start' )
[75]344             IF ( .NOT. humidity ) THEN
[97]345                IF ( ocean )  THEN
[388]346                   CALL diffusivities( prho, prho_reference )
[97]347                ELSE
348                   CALL diffusivities( pt, pt_reference )
349                ENDIF
[1]350             ELSE
[97]351                CALL diffusivities( vpt, pt_reference )
[1]352             ENDIF
353             CALL cpu_log( log_point(17), 'diffusivities', 'stop' )
[1015]354!
355!++          Statistics still require update of diffusivities on host
356             !$acc update host( kh, km )
[1]357
358          ENDIF
359
360       ENDDO   ! Intermediate step loop
361
362!
363!--    Increase simulation time and output times
[291]364       current_timestep_number    = current_timestep_number + 1
365       simulated_time             = simulated_time   + dt_3d
366       simulated_time_chr         = time_to_string( simulated_time )
367       time_since_reference_point = simulated_time - coupling_start_time
368
[1]369       IF ( simulated_time >= skip_time_data_output_av )  THEN
370          time_do_av         = time_do_av       + dt_3d
371       ENDIF
372       IF ( simulated_time >= skip_time_do2d_xy )  THEN
373          time_do2d_xy       = time_do2d_xy     + dt_3d
374       ENDIF
375       IF ( simulated_time >= skip_time_do2d_xz )  THEN
376          time_do2d_xz       = time_do2d_xz     + dt_3d
377       ENDIF
378       IF ( simulated_time >= skip_time_do2d_yz )  THEN
379          time_do2d_yz       = time_do2d_yz     + dt_3d
380       ENDIF
381       IF ( simulated_time >= skip_time_do3d    )  THEN
382          time_do3d          = time_do3d        + dt_3d
383       ENDIF
[410]384       DO  mid = 1, masks
385          IF ( simulated_time >= skip_time_domask(mid) )  THEN
386             time_domask(mid)= time_domask(mid) + dt_3d
387          ENDIF
388       ENDDO
[1]389       time_dvrp          = time_dvrp        + dt_3d
390       IF ( simulated_time >= skip_time_dosp )  THEN
391          time_dosp       = time_dosp        + dt_3d
392       ENDIF
393       time_dots          = time_dots        + dt_3d
[849]394       IF ( .NOT. first_call_lpm )  THEN
[1]395          time_dopts      = time_dopts       + dt_3d
396       ENDIF
397       IF ( simulated_time >= skip_time_dopr )  THEN
398          time_dopr       = time_dopr        + dt_3d
399       ENDIF
400       time_dopr_listing          = time_dopr_listing        + dt_3d
401       time_run_control   = time_run_control + dt_3d
402
403!
[102]404!--    Data exchange between coupled models
[291]405       IF ( coupling_mode /= 'uncoupled'  .AND.  run_coupled )  THEN
[102]406          time_coupling = time_coupling + dt_3d
[343]407
[108]408!
409!--       In case of model termination initiated by the local model
410!--       (terminate_coupled > 0), the coupler must be skipped because it would
411!--       cause an MPI intercomminucation hang.
412!--       If necessary, the coupler will be called at the beginning of the
413!--       next restart run.
414          DO WHILE ( time_coupling >= dt_coupling .AND. terminate_coupled == 0 )
[102]415             CALL surface_coupler
416             time_coupling = time_coupling - dt_coupling
417          ENDDO
418       ENDIF
419
420!
[46]421!--    Execute user-defined actions
422       CALL user_actions( 'after_integration' )
423
424!
[1]425!--    If Galilei transformation is used, determine the distance that the
426!--    model has moved so far
427       IF ( galilei_transformation )  THEN
428          advected_distance_x = advected_distance_x + u_gtrans * dt_3d
429          advected_distance_y = advected_distance_y + v_gtrans * dt_3d
430       ENDIF
431
432!
433!--    Check, if restart is necessary (because cpu-time is expiring or
434!--    because it is forced by user) and set stop flag
[108]435!--    This call is skipped if the remote model has already initiated a restart.
436       IF ( .NOT. terminate_run )  CALL check_for_restart
[1]437
438!
439!--    Carry out statistical analysis and output at the requested output times.
440!--    The MOD function is used for calculating the output time counters (like
441!--    time_dopr) in order to regard a possible decrease of the output time
442!--    interval in case of restart runs
443
444!
445!--    Set a flag indicating that so far no statistics have been created
446!--    for this time step
447       flow_statistics_called = .FALSE.
448
449!
450!--    If required, call flow_statistics for averaging in time
451       IF ( averaging_interval_pr /= 0.0  .AND.  &
452            ( dt_dopr - time_dopr ) <= averaging_interval_pr  .AND.  &
453            simulated_time >= skip_time_dopr )  THEN
454          time_dopr_av = time_dopr_av + dt_3d
455          IF ( time_dopr_av >= dt_averaging_input_pr )  THEN
456             do_sum = .TRUE.
457             time_dopr_av = MOD( time_dopr_av, &
458                                    MAX( dt_averaging_input_pr, dt_3d ) )
459          ENDIF
460       ENDIF
461       IF ( do_sum )  CALL flow_statistics
462
463!
[410]464!--    Sum-up 3d-arrays for later output of time-averaged 2d/3d/masked data
[1]465       IF ( averaging_interval /= 0.0  .AND.                                &
466            ( dt_data_output_av - time_do_av ) <= averaging_interval  .AND. &
467            simulated_time >= skip_time_data_output_av )                    &
468       THEN
469          time_do_sla = time_do_sla + dt_3d
470          IF ( time_do_sla >= dt_averaging_input )  THEN
471             CALL sum_up_3d_data
472             average_count_3d = average_count_3d + 1
473             time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) )
474          ENDIF
475       ENDIF
476
477!
478!--    Calculate spectra for time averaging
479       IF ( averaging_interval_sp /= 0.0  .AND.  &
480            ( dt_dosp - time_dosp ) <= averaging_interval_sp  .AND.  &
481            simulated_time >= skip_time_dosp )  THEN
482          time_dosp_av = time_dosp_av + dt_3d
483          IF ( time_dosp_av >= dt_averaging_input_pr )  THEN
484             CALL calc_spectra
485             time_dosp_av = MOD( time_dosp_av, &
486                                 MAX( dt_averaging_input_pr, dt_3d ) )
487          ENDIF
488       ENDIF
489
490!
491!--    Computation and output of run control parameters.
[1001]492!--    This is also done whenever perturbations have been imposed
[1]493       IF ( time_run_control >= dt_run_control  .OR.                     &
[1001]494            timestep_scheme(1:5) /= 'runge'  .OR.  disturbance_created ) &
[1]495       THEN
496          CALL run_control
497          IF ( time_run_control >= dt_run_control )  THEN
498             time_run_control = MOD( time_run_control, &
499                                     MAX( dt_run_control, dt_3d ) )
500          ENDIF
501       ENDIF
502
503!
504!--    Profile output (ASCII) on file
505       IF ( time_dopr_listing >= dt_dopr_listing )  THEN
506          CALL print_1d
507          time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, &
508                                                           dt_3d ) )
509       ENDIF
510
511!
512!--    Graphic output for PROFIL
513       IF ( time_dopr >= dt_dopr )  THEN
514          IF ( dopr_n /= 0 )  CALL data_output_profiles
515          time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) )
516          time_dopr_av = 0.0    ! due to averaging (see above)
517       ENDIF
518
519!
520!--    Graphic output for time series
521       IF ( time_dots >= dt_dots )  THEN
[48]522          CALL data_output_tseries
[1]523          time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) )
524       ENDIF
525
526!
527!--    Output of spectra (formatted for use with PROFIL), in case of no
528!--    time averaging, spectra has to be calculated before
529       IF ( time_dosp >= dt_dosp )  THEN
530          IF ( average_count_sp == 0 )  CALL calc_spectra
531          CALL data_output_spectra
532          time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) )
533       ENDIF
534
535!
536!--    2d-data output (cross-sections)
537       IF ( time_do2d_xy >= dt_do2d_xy )  THEN
538          CALL data_output_2d( 'xy', 0 )
539          time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) )
540       ENDIF
541       IF ( time_do2d_xz >= dt_do2d_xz )  THEN
542          CALL data_output_2d( 'xz', 0 )
543          time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) )
544       ENDIF
545       IF ( time_do2d_yz >= dt_do2d_yz )  THEN
546          CALL data_output_2d( 'yz', 0 )
547          time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) )
548       ENDIF
549
550!
551!--    3d-data output (volume data)
552       IF ( time_do3d >= dt_do3d )  THEN
553          CALL data_output_3d( 0 )
554          time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) )
555       ENDIF
556
557!
[410]558!--    masked data output
559       DO  mid = 1, masks
560          IF ( time_domask(mid) >= dt_domask(mid) )  THEN
561             CALL data_output_mask( 0 )
562             time_domask(mid) = MOD( time_domask(mid),  &
563                                     MAX( dt_domask(mid), dt_3d ) )
564          ENDIF
565       ENDDO
566
567!
568!--    Output of time-averaged 2d/3d/masked data
[1]569       IF ( time_do_av >= dt_data_output_av )  THEN
570          CALL average_3d_data
571          CALL data_output_2d( 'xy', 1 )
572          CALL data_output_2d( 'xz', 1 )
573          CALL data_output_2d( 'yz', 1 )
574          CALL data_output_3d( 1 )
[410]575          DO  mid = 1, masks
576             CALL data_output_mask( 1 )
577          ENDDO
[1]578          time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) )
579       ENDIF
580
581!
582!--    Output of particle time series
[253]583       IF ( particle_advection )  THEN
584          IF ( time_dopts >= dt_dopts  .OR. &
585               ( simulated_time >= particle_advection_start  .AND. &
[849]586                 first_call_lpm ) )  THEN
[253]587             CALL data_output_ptseries
588             time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
589          ENDIF
[1]590       ENDIF
591
592!
593!--    Output of dvrp-graphics (isosurface, particles, slicer)
594#if defined( __dvrp_graphics )
595       CALL DVRP_LOG_EVENT( -2, current_timestep_number-1 )
596#endif
597       IF ( time_dvrp >= dt_dvrp )  THEN
598          CALL data_output_dvrp
599          time_dvrp = MOD( time_dvrp, MAX( dt_dvrp, dt_3d ) )
600       ENDIF
601#if defined( __dvrp_graphics )
602       CALL DVRP_LOG_EVENT( 2, current_timestep_number )
603#endif
604
605!
606!--    If required, set the heat flux for the next time step at a random value
607       IF ( constant_heatflux  .AND.  random_heatflux )  CALL disturb_heatflux
608
609!
610!--    Execute user-defined actions
611       CALL user_actions( 'after_timestep' )
612
613       CALL cpu_log( log_point_s(10), 'timesteps', 'stop' )
614
[667]615
[1]616    ENDDO   ! time loop
617
618#if defined( __dvrp_graphics )
619    CALL DVRP_LOG_EVENT( -2, current_timestep_number )
620#endif
621
622 END SUBROUTINE time_integration
Note: See TracBrowser for help on using the repository browser.