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

Last change on this file since 404 was 398, checked in by raasch, 15 years ago

bugfix: exchange of ghost points for prho included

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