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

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

New:
---
ocean version including prognostic equation for salinity and equation of state for seawater. Routine buoyancy can be used with both temperature and density.
+ inipar-parameters bc_sa_t, bottom_salinityflux, ocean, sa_surface, sa_vertical_gradient, sa_vertical_gradient_level, top_salinityflux

advec_s_bc, average_3d_data, boundary_conds, buoyancy, check_parameters, data_output_2d, data_output_3d, diffusion_e, flow_statistics, header, init_grid, init_3d_model, modules, netcdf, parin, production_e, prognostic_equations, read_var_list, sum_up_3d_data, swap_timelevel, time_integration, user_interface, write_var_list, write_3d_binary

New:
eqn_state_seawater, init_ocean

Changed:


inipar-parameter use_pt_reference renamed use_reference

hydro_press renamed hyp, routine calc_mean_pt_profile renamed calc_mean_profile

format adjustments for the ocean version (run_control)

advec_particles, buoyancy, calc_liquid_water_content, check_parameters, diffusion_e, diffusivities, header, init_cloud_physics, modules, production_e, prognostic_equations, run_control

Errors:


Bugfix: height above topography instead of height above level k=0 is used for calculating the mixing length (diffusion_e and diffusivities).

Bugfix: error in boundary condition for TKE removed (advec_s_bc)

advec_s_bc, diffusion_e, prognostic_equations

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