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

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

further preliminary uncomplete changes for ocean version

  • Property svn:keywords set to Id
File size: 16.5 KB
RevLine 
[1]1 SUBROUTINE time_integration
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
[95]6! Ghostpoint exchange for salinity and density
[77]7!
8! Former revisions:
9! -----------------
10! $Id: time_integration.f90 95 2007-06-02 16:48:38Z raasch $
11!
[90]12! 87 2007-05-22 15:46:47Z raasch
13! var_hom renamed pr_palm
14!
[77]15! 75 2007-03-22 09:54:05Z raasch
[46]16! Move call of user_actions( 'after_integration' ) below increment of times
[63]17! and counters,
18! calls of prognostic_equations_.. changed to .._noopt, .._cache, and
[75]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
[1]22!
[3]23! RCS Log replace by Id keyword, revision history cleaned up
24!
[1]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.
[63]111          IF ( loop_optimization == 'vector' )  THEN
112             CALL prognostic_equations_vector
[1]113          ELSE
114             IF ( momentum_advec == 'ups-scheme'  .OR.  &
115                  scalar_advec == 'ups-scheme'   .OR.  &
116                  scalar_advec == 'bc-scheme' )        &
117             THEN
[63]118                CALL prognostic_equations_noopt
[1]119             ELSE
[63]120                CALL prognostic_equations_cache
[1]121             ENDIF
122          ENDIF
123
124!
125!--       Particle advection (only once during intermediate steps, because
126!--       it uses an Euler-step)
[63]127          IF ( particle_advection  .AND.                         &
128               simulated_time >= particle_advection_start  .AND. &
[1]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' )
[75]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 )
[95]152          IF ( ocean )  THEN
153             CALL exchange_horiz( sa_p )
154             CALL exchange_horiz( rho )
155          ENDIF
[75]156          IF ( humidity  .OR.  passive_scalar )  CALL exchange_horiz( q_p )
[1]157          IF ( cloud_droplets )  THEN
[75]158             CALL exchange_horiz( ql )
159             CALL exchange_horiz( ql_c )
160             CALL exchange_horiz( ql_v )
161             CALL exchange_horiz( ql_vp )
[1]162          ENDIF
163
164          CALL cpu_log( log_point(26), 'exchange-horiz-progn', 'stop' )
165
166!
167!--       Apply time filter in case of leap-frog timestep
168          IF ( tsc(2) == 2.0  .AND.  timestep_scheme(1:8) == 'leapfrog' )  THEN
169             CALL asselin_filter
170          ENDIF
171
172!
173!--       Boundary conditions for the prognostic quantities (except of the
174!--       velocities at the outflow in case of a non-cyclic lateral wall)
175          CALL boundary_conds( 'main' )
176
177!
[73]178!--       Swap the time levels in preparation for the next time step.
179          CALL swap_timelevel
180
181!
[1]182!--       Temperature offset must be imposed at cyclic boundaries in x-direction
183!--       when a sloping surface is used
184          IF ( sloping_surface )  THEN
185             IF ( nxl ==  0 )  pt(:,:,nxl-1) = pt(:,:,nxl-1) - pt_slope_offset
186             IF ( nxr == nx )  pt(:,:,nxr+1) = pt(:,:,nxr+1) + pt_slope_offset
187          ENDIF
188
189!
190!--       Impose a random perturbation on the horizontal velocity field
191          IF ( create_disturbances  .AND.  &
192               intermediate_timestep_count == intermediate_timestep_count_max )&
193          THEN
194             time_disturb = time_disturb + dt_3d
195             IF ( time_disturb >= dt_disturb )  THEN
[87]196                IF ( hom(nzb+5,1,pr_palm,0) < disturbance_energy_limit )  THEN
[75]197                   CALL disturb_field( nzb_u_inner, tend, u )
198                   CALL disturb_field( nzb_v_inner, tend, v )
[1]199                ELSEIF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
200!
201!--                Runs with a non-cyclic lateral wall need perturbations
202!--                near the inflow throughout the whole simulation
203                   dist_range = 1
[75]204                   CALL disturb_field( nzb_u_inner, tend, u )
205                   CALL disturb_field( nzb_v_inner, tend, v )
[1]206                   dist_range = 0
207                ENDIF
208                time_disturb = time_disturb - dt_disturb
209             ENDIF
210          ENDIF
211
212!
213!--       Reduce the velocity divergence via the equation for perturbation
214!--       pressure.
215          IF ( intermediate_timestep_count == intermediate_timestep_count_max &
216               .OR.  call_psolver_at_all_substeps )  THEN
217             CALL pres
218          ENDIF
219
220!
221!--       In case of a non-cyclic lateral wall, set the boundary conditions for
222!--       the velocities at the outflow
[75]223!          IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
224!             CALL boundary_conds( 'outflow_uvw' )
225!          ENDIF
[1]226
227!
228!--       If required, compute virtuell potential temperature
[75]229          IF ( humidity ) CALL compute_vpt
[1]230
231!
232!--       If required, compute liquid water content
233          IF ( cloud_physics ) CALL calc_liquid_water_content
234
235!
236!--       Compute the diffusion quantities
237          IF ( .NOT. constant_diffusion )  THEN
238
239!
240!--          First the vertical fluxes in the Prandtl layer are being computed
241             IF ( prandtl_layer )  THEN
242                CALL cpu_log( log_point(19), 'prandtl_fluxes', 'start' )
243                CALL prandtl_fluxes
244                CALL cpu_log( log_point(19), 'prandtl_fluxes', 'stop' )
245             ENDIF
246
247!
248!--          Compute the diffusion coefficients
249             CALL cpu_log( log_point(17), 'diffusivities', 'start' )
[75]250             IF ( .NOT. humidity ) THEN
[1]251                CALL diffusivities( pt )
252             ELSE
253                CALL diffusivities( vpt )
254             ENDIF
255             CALL cpu_log( log_point(17), 'diffusivities', 'stop' )
256
257          ENDIF
258
259       ENDDO   ! Intermediate step loop
260
261!
262!--    Increase simulation time and output times
263       current_timestep_number = current_timestep_number + 1
264       simulated_time     = simulated_time   + dt_3d
265       simulated_time_chr = time_to_string( simulated_time )
266       IF ( simulated_time >= skip_time_data_output_av )  THEN
267          time_do_av         = time_do_av       + dt_3d
268       ENDIF
269       IF ( simulated_time >= skip_time_do2d_xy )  THEN
270          time_do2d_xy       = time_do2d_xy     + dt_3d
271       ENDIF
272       IF ( simulated_time >= skip_time_do2d_xz )  THEN
273          time_do2d_xz       = time_do2d_xz     + dt_3d
274       ENDIF
275       IF ( simulated_time >= skip_time_do2d_yz )  THEN
276          time_do2d_yz       = time_do2d_yz     + dt_3d
277       ENDIF
278       IF ( simulated_time >= skip_time_do3d    )  THEN
279          time_do3d          = time_do3d        + dt_3d
280       ENDIF
281       time_dvrp          = time_dvrp        + dt_3d
282       IF ( simulated_time >= skip_time_dosp )  THEN
283          time_dosp       = time_dosp        + dt_3d
284       ENDIF
285       time_dots          = time_dots        + dt_3d
286       IF ( .NOT. first_call_advec_particles )  THEN
287          time_dopts      = time_dopts       + dt_3d
288       ENDIF
289       IF ( simulated_time >= skip_time_dopr )  THEN
290          time_dopr       = time_dopr        + dt_3d
291       ENDIF
292       time_dopr_listing          = time_dopr_listing        + dt_3d
293       time_run_control   = time_run_control + dt_3d
294
295!
[46]296!--    Execute user-defined actions
297       CALL user_actions( 'after_integration' )
298
299!
[1]300!--    If Galilei transformation is used, determine the distance that the
301!--    model has moved so far
302       IF ( galilei_transformation )  THEN
303          advected_distance_x = advected_distance_x + u_gtrans * dt_3d
304          advected_distance_y = advected_distance_y + v_gtrans * dt_3d
305       ENDIF
306
307!
308!--    Check, if restart is necessary (because cpu-time is expiring or
309!--    because it is forced by user) and set stop flag
310       CALL check_for_restart
311
312!
313!--    Carry out statistical analysis and output at the requested output times.
314!--    The MOD function is used for calculating the output time counters (like
315!--    time_dopr) in order to regard a possible decrease of the output time
316!--    interval in case of restart runs
317
318!
319!--    Set a flag indicating that so far no statistics have been created
320!--    for this time step
321       flow_statistics_called = .FALSE.
322
323!
324!--    If required, call flow_statistics for averaging in time
325       IF ( averaging_interval_pr /= 0.0  .AND.  &
326            ( dt_dopr - time_dopr ) <= averaging_interval_pr  .AND.  &
327            simulated_time >= skip_time_dopr )  THEN
328          time_dopr_av = time_dopr_av + dt_3d
329          IF ( time_dopr_av >= dt_averaging_input_pr )  THEN
330             do_sum = .TRUE.
331             time_dopr_av = MOD( time_dopr_av, &
332                                    MAX( dt_averaging_input_pr, dt_3d ) )
333          ENDIF
334       ENDIF
335       IF ( do_sum )  CALL flow_statistics
336
337!
338!--    Sum-up 3d-arrays for later output of time-averaged data
339       IF ( averaging_interval /= 0.0  .AND.                                &
340            ( dt_data_output_av - time_do_av ) <= averaging_interval  .AND. &
341            simulated_time >= skip_time_data_output_av )                    &
342       THEN
343          time_do_sla = time_do_sla + dt_3d
344          IF ( time_do_sla >= dt_averaging_input )  THEN
345             CALL sum_up_3d_data
346             average_count_3d = average_count_3d + 1
347             time_do_sla = MOD( time_do_sla, MAX( dt_averaging_input, dt_3d ) )
348          ENDIF
349       ENDIF
350
351!
352!--    Calculate spectra for time averaging
353       IF ( averaging_interval_sp /= 0.0  .AND.  &
354            ( dt_dosp - time_dosp ) <= averaging_interval_sp  .AND.  &
355            simulated_time >= skip_time_dosp )  THEN
356          time_dosp_av = time_dosp_av + dt_3d
357          IF ( time_dosp_av >= dt_averaging_input_pr )  THEN
358             CALL calc_spectra
359             time_dosp_av = MOD( time_dosp_av, &
360                                 MAX( dt_averaging_input_pr, dt_3d ) )
361          ENDIF
362       ENDIF
363
364!
365!--    Computation and output of run control parameters.
366!--    This is also done whenever the time step has changed or perturbations
367!--    have been imposed
368       IF ( time_run_control >= dt_run_control  .OR.                     &
369            ( dt_changed  .AND.  timestep_scheme(1:5) /= 'runge' )  .OR. &
370            disturbance_created )                                        &
371       THEN
372          CALL run_control
373          IF ( time_run_control >= dt_run_control )  THEN
374             time_run_control = MOD( time_run_control, &
375                                     MAX( dt_run_control, dt_3d ) )
376          ENDIF
377       ENDIF
378
379!
380!--    Profile output (ASCII) on file
381       IF ( time_dopr_listing >= dt_dopr_listing )  THEN
382          CALL print_1d
383          time_dopr_listing = MOD( time_dopr_listing, MAX( dt_dopr_listing, &
384                                                           dt_3d ) )
385       ENDIF
386
387!
388!--    Graphic output for PROFIL
389       IF ( time_dopr >= dt_dopr )  THEN
390          IF ( dopr_n /= 0 )  CALL data_output_profiles
391          time_dopr = MOD( time_dopr, MAX( dt_dopr, dt_3d ) )
392          time_dopr_av = 0.0    ! due to averaging (see above)
393       ENDIF
394
395!
396!--    Graphic output for time series
397       IF ( time_dots >= dt_dots )  THEN
[48]398          CALL data_output_tseries
[1]399          time_dots = MOD( time_dots, MAX( dt_dots, dt_3d ) )
400       ENDIF
401
402!
403!--    Output of spectra (formatted for use with PROFIL), in case of no
404!--    time averaging, spectra has to be calculated before
405       IF ( time_dosp >= dt_dosp )  THEN
406          IF ( average_count_sp == 0 )  CALL calc_spectra
407          CALL data_output_spectra
408          time_dosp = MOD( time_dosp, MAX( dt_dosp, dt_3d ) )
409       ENDIF
410
411!
412!--    2d-data output (cross-sections)
413       IF ( time_do2d_xy >= dt_do2d_xy )  THEN
414          CALL data_output_2d( 'xy', 0 )
415          time_do2d_xy = MOD( time_do2d_xy, MAX( dt_do2d_xy, dt_3d ) )
416       ENDIF
417       IF ( time_do2d_xz >= dt_do2d_xz )  THEN
418          CALL data_output_2d( 'xz', 0 )
419          time_do2d_xz = MOD( time_do2d_xz, MAX( dt_do2d_xz, dt_3d ) )
420       ENDIF
421       IF ( time_do2d_yz >= dt_do2d_yz )  THEN
422          CALL data_output_2d( 'yz', 0 )
423          time_do2d_yz = MOD( time_do2d_yz, MAX( dt_do2d_yz, dt_3d ) )
424       ENDIF
425
426!
427!--    3d-data output (volume data)
428       IF ( time_do3d >= dt_do3d )  THEN
429          CALL data_output_3d( 0 )
430          time_do3d = MOD( time_do3d, MAX( dt_do3d, dt_3d ) )
431       ENDIF
432
433!
434!--    Output of time-averaged 2d/3d-data
435       IF ( time_do_av >= dt_data_output_av )  THEN
436          CALL average_3d_data
437          CALL data_output_2d( 'xy', 1 )
438          CALL data_output_2d( 'xz', 1 )
439          CALL data_output_2d( 'yz', 1 )
440          CALL data_output_3d( 1 )
441          time_do_av = MOD( time_do_av, MAX( dt_data_output_av, dt_3d ) )
442       ENDIF
443
444!
445!--    Output of particle time series
446       IF ( time_dopts >= dt_dopts  .OR. &
447            ( simulated_time >= particle_advection_start  .AND. &
448              first_call_advec_particles ) )  THEN
449          CALL data_output_ptseries
450          time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
451       ENDIF
452
453!
454!--    Output of dvrp-graphics (isosurface, particles, slicer)
455#if defined( __dvrp_graphics )
456       CALL DVRP_LOG_EVENT( -2, current_timestep_number-1 )
457#endif
458       IF ( time_dvrp >= dt_dvrp )  THEN
459          CALL data_output_dvrp
460          time_dvrp = MOD( time_dvrp, MAX( dt_dvrp, dt_3d ) )
461       ENDIF
462#if defined( __dvrp_graphics )
463       CALL DVRP_LOG_EVENT( 2, current_timestep_number )
464#endif
465
466!
467!--    If required, set the heat flux for the next time step at a random value
468       IF ( constant_heatflux  .AND.  random_heatflux )  CALL disturb_heatflux
469
470!
471!--    Execute user-defined actions
472       CALL user_actions( 'after_timestep' )
473
474       CALL cpu_log( log_point_s(10), 'timesteps', 'stop' )
475
476    ENDDO   ! time loop
477
478#if defined( __dvrp_graphics )
479    CALL DVRP_LOG_EVENT( -2, current_timestep_number )
480#endif
481
482 END SUBROUTINE time_integration
Note: See TracBrowser for help on using the repository browser.