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

Last change on this file since 332 was 291, checked in by raasch, 16 years ago

changes for coupling with independent precursor runs; z_i calculation with Sullivan criterion

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