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

Last change on this file since 597 was 484, checked in by raasch, 15 years ago

typo in file headers removed

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