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

Last change on this file since 635 was 623, checked in by raasch, 14 years ago

last commit documented

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