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

Last change on this file since 92 was 90, checked in by raasch, 17 years ago

New:
---
Calculation and output of user-defined profiles. New &userpar parameters data_output_pr_user and max_pr_user.

check_parameters, flow_statistics, modules, parin, read_var_list, user_interface, write_var_list

Changed:


Division through dt_3d replaced by multiplication of the inverse. For performance optimisation, this is done in the loop calculating the divergence instead of using a seperate loop. (pres.f90) var_hom and var_sum renamed pr_palm.

data_output_profiles, flow_statistics, init_3d_model, modules, parin, pres, read_var_list, run_control, time_integration

Errors:


Bugfix: work_fft*_vec removed from some PRIVATE-declarations (poisfft).

Bugfix: field_chr renamed field_char (user_interface).

Bugfix: output of use_upstream_for_tke (header).

header, poisfft, user_interface

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