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

Last change on this file since 790 was 790, checked in by raasch, 12 years ago

Bugfix for output of mean particle radius + preliminary works for implementing the Wang collision kernel

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