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

Last change on this file since 100 was 98, checked in by raasch, 17 years ago

updating comments and rc-file

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