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

Last change on this file since 803 was 791, checked in by raasch, 13 years ago

last commit documented

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