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

Last change on this file since 359 was 348, checked in by maronga, 15 years ago

bugfix for coupled restart runs

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