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

Last change on this file since 260 was 253, checked in by raasch, 16 years ago

small adjustments for NEC at RIAM; bugfix concerning output of particle timeseries

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