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

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

last commit documented

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