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

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

preliminary update of bugfixes and extensions for non-cyclic BCs

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