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

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

Changed:


Original routine advec_particles split into several new subroutines and renamed
lpm.
init_particles renamed lpm_init
user_advec_particles renamed user_lpm_advec,
particle_boundary_conds renamed lpm_boundary_conds,
set_particle_attributes renamed lpm_set_attributes,
user_init_particles renamed user_lpm_init,
user_particle_attributes renamed user_lpm_set_attributes
(Makefile, lpm_droplet_collision, lpm_droplet_condensation, init_3d_model, modules, palm, read_var_list, time_integration, write_var_list, deleted: advec_particles, init_particles, particle_boundary_conds, set_particle_attributes, user_advec_particles, user_init_particles, user_particle_attributes, new: lpm, lpm_advec, lpm_boundary_conds, lpm_calc_liquid_water_content, lpm_data_output_particles, lpm_droplet_collision, lpm_drollet_condensation, lpm_exchange_horiz, lpm_extend_particle_array, lpm_extend_tails, lpm_extend_tail_array, lpm_init, lpm_init_sgs_tke, lpm_pack_arrays, lpm_read_restart_file, lpm_release_set, lpm_set_attributes, lpm_sort_arrays, lpm_write_exchange_statistics, lpm_write_restart_file, user_lpm_advec, user_lpm_init, user_lpm_set_attributes

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