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

Last change on this file since 1060 was 1054, checked in by hoffmann, 12 years ago

last commit documented

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