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

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

New:
---
Allows runs for a coupled atmosphere-ocean LES,
coupling frequency is controlled by new d3par-parameter dt_coupling,
the coupling mode (atmosphere_to_ocean or ocean_to_atmosphere) for the
respective processes is read from environment variable coupling_mode,
which is set by the mpiexec-command,
communication between the two models is done using the intercommunicator
comm_inter,
local files opened by the ocean model get the additional suffic "_O".
Assume saturation at k=nzb_s_inner(j,i) for atmosphere coupled to ocean.

A momentum flux can be set as top boundary condition using the new
inipar parameter top_momentumflux_u|v.

Non-cyclic boundary conditions can be used along all horizontal directions.

Quantities w*p* and w"e can be output as vertical profiles.

Initial profiles are reset to constant profiles in case that initializing_actions /= 'set_constant_profiles'. (init_rankine)

Optionally calculate km and kh from initial TKE e_init.

Changed:


Remaining variables iran changed to iran_part (advec_particles, init_particles).

In case that the presure solver is not called for every Runge-Kutta substep
(call_psolver_at_all_substeps = .F.), it is called after the first substep
instead of the last. In that case, random perturbations are also added to the
velocity field after the first substep.

Initialization of km,kh = 0.00001 for ocean = .T. (for ocean = .F. it remains 0.01).

Allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of cloud_physics = .T.).

Errors:


Bugs from code parts for non-cyclic boundary conditions are removed: loops for
u and v are starting from index nxlu, nysv, respectively. The radiation boundary
condition is used for every Runge-Kutta substep. Velocity phase speeds for
the radiation boundary conditions are calculated for the first Runge-Kutta
substep only and reused for the further substeps. New arrays c_u, c_v, and c_w
are defined for this purpose. Several index errors are removed from the
radiation boundary condition code parts. Upper bounds for calculating
u_0 and v_0 (in production_e) are nxr+1 and nyn+1 because otherwise these
values are not available in case of non-cyclic boundary conditions.

+dots_num_palm in module user, +module netcdf_control in user_init (both in user_interface)

Bugfix: wrong sign removed from the buoyancy production term in the case use_reference = .T. (production_e)

Bugfix: Error message concerning output of particle concentration (pc) modified (check_parameters).

Bugfix: Rayleigh damping for ocean fixed.

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