source: palm/trunk/SOURCE/timestep.f90 @ 753

Last change on this file since 753 was 708, checked in by raasch, 14 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 15.9 KB
RevLine 
[1]1 SUBROUTINE timestep
2
3!------------------------------------------------------------------------------!
[258]4! Current revisions:
[1]5! -----------------
[708]6!
[316]7!
[1]8! Former revisions:
9! -----------------
[3]10! $Id: timestep.f90 708 2011-03-29 12:34:54Z raasch $
[110]11!
[708]12! 707 2011-03-29 11:39:40Z raasch
13! bc_lr/ns replaced by bc_lr/ns_cyc
14!
[668]15! 667 2010-12-23 12:06:00Z suehring/gryschka
16! Exchange of terminate_coupled between ocean and atmosphere via PE0
17! Minimum grid spacing dxyz2_min(k) is now calculated using dzw instead of dzu
18!
[623]19! 622 2010-12-10 08:08:13Z raasch
20! optional barriers included in order to speed up collective operations
21!
[392]22! 343 2009-06-24 12:59:09Z maronga
23! Additional timestep criterion in case of simulations with plant canopy
24! Output of messages replaced by message handling routine.
25!
[226]26! 222 2009-01-12 16:04:16Z letzel
27! Implementation of a MPI-1 Coupling: replaced myid with target_id
28! Bugfix for nonparallel execution
29!
[110]30! 108 2007-08-24 15:10:38Z letzel
31! modifications to terminate coupled runs
32!
[3]33! RCS Log replace by Id keyword, revision history cleaned up
34!
[1]35! Revision 1.21  2006/02/23 12:59:44  raasch
36! nt_anz renamed current_timestep_number
37!
38! Revision 1.1  1997/08/11 06:26:19  raasch
39! Initial revision
40!
41!
42! Description:
43! ------------
44! Compute the time step under consideration of the FCL and diffusion criterion.
45!------------------------------------------------------------------------------!
46
47    USE arrays_3d
48    USE control_parameters
49    USE cpulog
50    USE grid_variables
51    USE indices
52    USE interfaces
53    USE pegrid
54    USE statistics
55
56    IMPLICIT NONE
57
58    INTEGER ::  i, j, k
59
[318]60    REAL ::  div, dt_diff, dt_diff_l, dt_plant_canopy,                 &
61             dt_plant_canopy_l,                                        &
62             dt_plant_canopy_u, dt_plant_canopy_v, dt_plant_canopy_w,  & 
63             dt_u, dt_v, dt_w, lad_max, percent_change,                &
[316]64             u_gtrans_l, vabs_max, value, v_gtrans_l
[1]65
66    REAL, DIMENSION(2)         ::  uv_gtrans, uv_gtrans_l
67    REAL, DIMENSION(nzb+1:nzt) ::  dxyz2_min
68
[667]69
70
[1]71    CALL cpu_log( log_point(12), 'calculate_timestep', 'start' )
72
73!
74!-- Determine the maxima of the velocity components.
[667]75    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, u, 'abs', &
[1]76                         u_max, u_max_ijk )
[667]77    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, v, 'abs', &
[1]78                         v_max, v_max_ijk )
[667]79    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, w, 'abs', &
[1]80                         w_max, w_max_ijk )
81
82!
83!-- In case maxima of the horizontal velocity components have been found at the
84!-- bottom boundary (k=nzb), the corresponding maximum at level k=1 is chosen
85!-- if the Dirichlet-boundary condition ('mirror') has been set. This is
86!-- necessary, because otherwise in case of Galilei-transform a far too large
87!-- velocity (having the respective opposite sign) would be used for the time
88!-- step determination (almost double the mean flow velocity).
89    IF ( ibc_uv_b == 0 )  THEN
90       IF ( u_max_ijk(1) == nzb )  THEN
91          u_max        = -u_max
92          u_max_ijk(1) = nzb + 1
93       ENDIF
94       IF ( v_max_ijk(1) == nzb )  THEN
95          v_max        = -v_max
96          v_max_ijk(1) = nzb + 1
97       ENDIF
98    ENDIF
99
100!
101!-- In case of Galilei-transform not using the geostrophic wind as translation
102!-- velocity, compute the volume-averaged horizontal velocity components, which
103!-- will then be subtracted from the horizontal wind for the time step and
104!-- horizontal advection routines.
105    IF ( galilei_transformation  .AND. .NOT.  use_ug_for_galilei_tr )  THEN
106       IF ( flow_statistics_called )  THEN
107!
108!--       Horizontal averages already existent, just need to average them
109!--       vertically.
110          u_gtrans = 0.0
111          v_gtrans = 0.0
112          DO  k = nzb+1, nzt
113             u_gtrans = u_gtrans + hom(k,1,1,0)
114             v_gtrans = v_gtrans + hom(k,1,2,0)
115          ENDDO
116          u_gtrans = u_gtrans / REAL( nzt - nzb )
117          v_gtrans = v_gtrans / REAL( nzt - nzb )
118       ELSE
119!
120!--       Averaging over the entire model domain.
121          uv_gtrans_l = 0.0
122          DO  i = nxl, nxr
123             DO  j = nys, nyn
124                DO  k = nzb+1, nzt
125                   uv_gtrans_l(1) = uv_gtrans_l(1) + u(k,j,i)
126                   uv_gtrans_l(2) = uv_gtrans_l(2) + v(k,j,i)
127                ENDDO
128             ENDDO
129          ENDDO
130          uv_gtrans_l = uv_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) )
131#if defined( __parallel )
[622]132          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[1]133          CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM, &
134                              comm2d, ierr )
135          u_gtrans = uv_gtrans(1) / REAL( numprocs )
136          v_gtrans = uv_gtrans(2) / REAL( numprocs )
137#else
138          u_gtrans = uv_gtrans_l(1)
139          v_gtrans = uv_gtrans_l(2)
140#endif
141       ENDIF
142    ENDIF
143
144    IF ( .NOT. dt_fixed )  THEN
145!
146!--    Variable time step:
147!
148!--    For each component, compute the maximum time step according to the
149!--    FCL-criterion.
150       dt_u = dx / ( ABS( u_max - u_gtrans ) + 1.0E-10 )
151       dt_v = dy / ( ABS( v_max - v_gtrans ) + 1.0E-10 )
152       dt_w = dzu(MAX( 1, w_max_ijk(1) )) / ( ABS( w_max ) + 1.0E-10 )
153
154!
155!--    Compute time step according to the diffusion criterion.
156!--    First calculate minimum grid spacing which only depends on index k
157!--    Note: also at k=nzb+1 a full grid length is being assumed, although
158!--          in the Prandtl-layer friction term only dz/2 is used.
159!--          Experience from the old model seems to justify this.
160       dt_diff_l = 999999.0
161
162       DO  k = nzb+1, nzt
[667]163           dxyz2_min(k) = MIN( dx2, dy2, dzw(k)*dzw(k) ) * 0.125
[1]164       ENDDO
165
166!$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l)
167!$OMP DO
168       DO  i = nxl, nxr
169          DO  j = nys, nyn
170             DO  k = nzb+1, nzt
171                value = dxyz2_min(k) / ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 )
172
173                dt_diff_l = MIN( value, dt_diff_l )
174             ENDDO
175          ENDDO
176       ENDDO
177!$OMP END PARALLEL
178#if defined( __parallel )
[622]179       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[1]180       CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, &
181                           ierr )
182#else
183       dt_diff = dt_diff_l
184#endif
185
186!
187!--    In case of non-cyclic lateral boundaries, the diffusion time step
188!--    may be further restricted by the lateral damping layer (damping only
189!--    along x and y)
[707]190       IF ( .NOT. bc_lr_cyc )  THEN
[1]191          dt_diff = MIN( dt_diff, 0.125 * dx2 / ( km_damp_max + 1E-20 ) )
[707]192       ELSEIF ( .NOT. bc_ns_cyc )  THEN
[1]193          dt_diff = MIN( dt_diff, 0.125 * dy2 / ( km_damp_max + 1E-20 ) )
194       ENDIF
195
196!
[316]197!--    Additional timestep criterion with plant canopies:
198!--    it is not allowed to extract more than the available momentum
199       IF ( plant_canopy ) THEN
[318]200
201          dt_plant_canopy_l = 0.0
202          DO  i = nxl, nxr
203             DO  j = nys, nyn
204                DO k = nzb+1, nzt
205                   dt_plant_canopy_u = cdc(k,j,i) * lad_u(k,j,i) *  &
206                                       SQRT(     u(k,j,i)**2     +  &
207                                             ( ( v(k,j,i-1)      +  &
208                                                 v(k,j,i)        +  &
209                                                 v(k,j+1,i)      +  &
210                                                 v(k,j+1,i-1) )     &
211                                               / 4.0 )**2        +  &
212                                             ( ( w(k-1,j,i-1)    +  &
213                                                 w(k-1,j,i)      +  &
214                                                 w(k,j,i-1)      +  &
215                                                 w(k,j,i) )         &
216                                                 / 4.0 )**2 ) 
217                   IF ( dt_plant_canopy_u > dt_plant_canopy_l ) THEN
218                      dt_plant_canopy_l = dt_plant_canopy_u 
219                   ENDIF
220                   dt_plant_canopy_v = cdc(k,j,i) * lad_v(k,j,i) *  &
221                                       SQRT( ( ( u(k,j-1,i)      +  &
222                                                 u(k,j-1,i+1)    +  &
223                                                 u(k,j,i)        +  &
224                                                 u(k,j,i+1) )       &
225                                               / 4.0 )**2        +  &
226                                                 v(k,j,i)**2     +  &
227                                             ( ( w(k-1,j-1,i)    +  &
228                                                 w(k-1,j,i)      +  &
229                                                 w(k,j-1,i)      +  &
230                                                 w(k,j,i) )         &
231                                                 / 4.0 )**2 ) 
232                   IF ( dt_plant_canopy_v > dt_plant_canopy_l ) THEN
233                      dt_plant_canopy_l = dt_plant_canopy_v
234                   ENDIF                   
235                   dt_plant_canopy_w = cdc(k,j,i) * lad_w(k,j,i) *  &
236                                       SQRT( ( ( u(k,j,i)        +  &
237                                                 u(k,j,i+1)      +  &
238                                                 u(k+1,j,i)      +  &
239                                                 u(k+1,j,i+1) )     &
240                                               / 4.0 )**2        +  &
241                                             ( ( v(k,j,i)        +  &
242                                                 v(k,j+1,i)      +  &
243                                                 v(k+1,j,i)      +  &
244                                                 v(k+1,j+1,i) )     &
245                                               / 4.0 )**2        +  &
246                                                 w(k,j,i)**2 )     
247                   IF ( dt_plant_canopy_w > dt_plant_canopy_l ) THEN
248                      dt_plant_canopy_l = dt_plant_canopy_w
249                   ENDIF
250                ENDDO
251             ENDDO
252          ENDDO 
253
254          IF ( dt_plant_canopy_l > 0.0 ) THEN
[320]255!
256!--          Invert dt_plant_canopy_l and apply a security timestep factor 0.1
[318]257             dt_plant_canopy_l = 0.1 / dt_plant_canopy_l
[320]258          ELSE
259!
260!--          In case of inhomogeneous plant canopy, some processors may have no
261!--          canopy at all. Then use dt_max as dummy instead.
262             dt_plant_canopy_l = dt_max
[318]263          ENDIF
[320]264
[316]265!
[318]266!--       Determine the global minumum
267#if defined( __parallel )
[622]268          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
[318]269          CALL MPI_ALLREDUCE( dt_plant_canopy_l, dt_plant_canopy, 1, MPI_REAL,  &
270                              MPI_MIN, comm2d, ierr )
271#else
272          dt_plant_canopy = dt_plant_canopy_l
273#endif
[316]274
275       ELSE
276!
277!--       Use dt_diff as dummy value to avoid extra IF branches further below
278          dt_plant_canopy = dt_diff
279
280       ENDIF
281
282!
283!--    The time step is the minimum of the 3-4 components and the diffusion time
[1]284!--    step minus a reduction to be on the safe side. Factor 0.5 is necessary
285!--    since the leap-frog scheme always progresses by 2 * delta t.
286!--    The user has to set the cfl_factor small enough to ensure that the
287!--    divergences do not become too large.
288!--    The time step must not exceed the maximum allowed value.
289       IF ( timestep_scheme(1:5) == 'runge' )  THEN
[316]290          dt_3d = cfl_factor * MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w )
[1]291       ELSE
[316]292          dt_3d = cfl_factor * 0.5 *  &
293                  MIN( dt_diff, dt_plant_canopy, dt_u, dt_v, dt_w )
[1]294       ENDIF
295       dt_3d = MIN( dt_3d, dt_max )
296
297!
298!--    Remember the restricting time step criterion for later output.
[316]299       IF ( MIN( dt_u, dt_v, dt_w ) < MIN( dt_diff, dt_plant_canopy ) )  THEN
[1]300          timestep_reason = 'A'
[316]301       ELSEIF ( dt_plant_canopy < dt_diff )  THEN
302          timestep_reason = 'C'
[1]303       ELSE
304          timestep_reason = 'D'
305       ENDIF
306
307!
308!--    Set flag if the time step becomes too small.
309       IF ( dt_3d < ( 0.00001 * dt_max ) )  THEN
310          stop_dt = .TRUE.
[108]311
[320]312          WRITE( message_string, * ) 'Time step has reached minimum limit.',  &
313               '&dt              = ', dt_3d, ' s  Simulation is terminated.', &
314               '&old_dt          = ', old_dt, ' s',                           &
315               '&dt_u            = ', dt_u, ' s',                             &
316               '&dt_v            = ', dt_v, ' s',                             &
317               '&dt_w            = ', dt_w, ' s',                             &
318               '&dt_diff         = ', dt_diff, ' s',                          &
319               '&dt_plant_canopy = ', dt_plant_canopy, ' s',                  &
320               '&u_max   = ', u_max, ' m/s   k=', u_max_ijk(1),               &
321               '  j=', u_max_ijk(2), '  i=', u_max_ijk(3),                    &
322               '&v_max   = ', v_max, ' m/s   k=', v_max_ijk(1),               &
323               '  j=', v_max_ijk(2), '  i=', v_max_ijk(3),                    &
324               '&w_max   = ', w_max, ' m/s   k=', w_max_ijk(1),               &
325               '  j=', w_max_ijk(2), '  i=', w_max_ijk(3)
[258]326          CALL message( 'timestep', 'PA0312', 0, 1, 0, 6, 0 )
[108]327!
328!--       In case of coupled runs inform the remote model of the termination
329!--       and its reason, provided the remote model has not already been
330!--       informed of another termination reason (terminate_coupled > 0) before.
[222]331#if defined( __parallel )
[108]332          IF ( coupling_mode /= 'uncoupled' .AND. terminate_coupled == 0 )  THEN
333             terminate_coupled = 2
[667]334             IF ( myid == 0 ) THEN
335                CALL MPI_SENDRECV( &
336                     terminate_coupled,        1, MPI_INTEGER, target_id,  0, &
337                     terminate_coupled_remote, 1, MPI_INTEGER, target_id,  0, &
338                     comm_inter, status, ierr )
339             ENDIF
340             CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0, comm2d, ierr)
[108]341          ENDIF
[222]342#endif
[1]343       ENDIF
344
345!
346!--    Ensure a smooth value (two significant digits) of the timestep. For
347!--    other schemes than Runge-Kutta, the following restrictions appear:
348!--    The current timestep is only then changed, if the change relative to
349!--    its previous value exceeds +5 % or -2 %. In case of a timestep
350!--    reduction, at least 30 iterations have to be performed before a timestep
351!--    enlargement is permitted again.
352       percent_change = dt_3d / old_dt - 1.0
353       IF ( percent_change > 0.05  .OR.  percent_change < -0.02  .OR. &
354            timestep_scheme(1:5) == 'runge' )  THEN
355
356!
357!--       Time step enlargement by no more than 2 %.
358          IF ( percent_change > 0.0  .AND.  simulated_time /= 0.0  .AND. &
359               timestep_scheme(1:5) /= 'runge' )  THEN
360             dt_3d = 1.02 * old_dt
361          ENDIF
362
363!
364!--       A relatively smooth value of the time step is ensured by taking
365!--       only the first two significant digits.
366          div = 1000.0
367          DO  WHILE ( dt_3d < div )
368             div = div / 10.0
369          ENDDO
370          dt_3d = NINT( dt_3d * 100.0 / div ) * div / 100.0
371
372!
373!--       Now the time step can be adjusted.
374          IF ( percent_change < 0.0  .OR.  timestep_scheme(1:5) == 'runge' ) &
375          THEN
376!
377!--          Time step reduction.
378             old_dt = dt_3d
379             dt_changed = .TRUE.
380          ELSE
381!
382!--          For other timestep schemes , the time step is only enlarged
383!--          after at least 30 iterations since the previous time step
384!--          change or, of course, after model initialization.
385             IF ( current_timestep_number >= last_dt_change + 30  .OR. &
386                  simulated_time == 0.0 )  THEN
387                old_dt = dt_3d
388                dt_changed = .TRUE.
389             ELSE
390                dt_3d = old_dt
391                dt_changed = .FALSE.
392             ENDIF
393
394          ENDIF
395       ELSE
396!
397!--       No time step change since the difference is too small.
398          dt_3d = old_dt
399          dt_changed = .FALSE.
400       ENDIF
401
402       IF ( dt_changed )  last_dt_change = current_timestep_number
403
404    ENDIF
405    CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' )
406
407 END SUBROUTINE timestep
Note: See TracBrowser for help on using the repository browser.