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

Last change on this file since 1 was 1, checked in by raasch, 17 years ago

Initial repository layout and content

File size: 11.8 KB
Line 
1 SUBROUTINE timestep
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Log: timestep.f90,v $
11! Revision 1.21  2006/02/23 12:59:44  raasch
12! nt_anz renamed current_timestep_number
13!
14! Revision 1.20  2005/04/23 09:45:38  raasch
15! fcl_factor renamed cfl_factor
16!
17! Revision 1.19  2005/03/26 21:11:53  raasch
18! In case of non-cyclic lateral boundary conditions, the increased diffusivity
19! in the lateral damping layer is regarded in the minimum timestep calculation.
20!
21! Revision 1.18  2004/01/30 10:40:02  raasch
22! Timestep increment limitation in case of Runge-Kutta removed
23!
24! Revision 1.17  2004/01/28 15:34:03  raasch
25! Setting of steering factors for the prognostic equations moved to new
26! routine timestep_scheme_steering now called in routine leap_frog.
27!
28! Revision 1.16  2003/03/16 09:49:47  raasch
29! Two underscores (_) are placed in front of all define-strings
30!
31! Revision 1.15  2002/12/19 16:19:54  raasch
32! Two mpi_allreduce calls for determining mean horizontal velocity reduced to
33! one call. Calculation of minimum timestep for diffusion optimized by
34! determining minimum grid spacing in a single k loop.
35!
36! Revision 1.14  2001/03/30 07:55:12  raasch
37! Translation of remaining German identifiers (variables, subroutines, etc.)
38!
39! Revision 1.13  2001/01/30 22:10:11  letzel
40! All comments translated into English.
41!
42! Revision 1.12  2001/01/22 08:17:33  raasch
43! Module test_variables removed
44!
45! Revision 1.11  2000/01/10 10:10:30  raasch
46! Translationsgeschwindigkeit fuer Galilei-Transformation wird nur noch
47! bedingt berechnet
48!
49! Revision 1.10  1999/02/05 09:19:06  raasch
50! Zeitschrittsteuerung beruecksichtigt jetzt auch Wahl des reinen Euler-
51! Verfahrens
52!
53! Revision 1.9  1998/09/22 17:30:41  raasch
54! Ausfuehrlichere Informationen bei Abbruch durch Unterschreitung des
55! minimalen Zeitschritts
56!
57! Revision 1.8  1998/07/06 12:33:52  raasch
58! Korrektur der ermittelten horizontalen Geschwindigkeitsmaxima, falls diese
59! bei k=0 gefunden wurden (wegen Spiegelungsrandbedingung treten dort negative
60! Werte auf, die bei Galilei-Transformation den Zeitschritt erheblich
61! einschraenken). + USE test_variables
62!
63! Revision 1.7  1998/04/21 15:57:46  raasch
64! Implementierung der Galilei-Transformation
65!
66! Revision 1.6  1998/03/25 13:57:20  raasch
67! dt in dt_3d umbenannt
68!
69! Revision 1.5  1998/03/09 16:25:10  raasch
70! Schoenheitskorrekturen
71!
72! Revision 1.4  1998/03/03 08:02:15  raasch
73! Leap-Frog auch nach Zeitschrittaenderung erlaubt
74!
75! Revision 1.3  1997/09/12 06:30:27  raasch
76! Einbau des Diffusions-Kriteriums
77!
78! Revision 1.2  1997/08/26 06:34:55  raasch
79! Fehler bei w_max_ijk - Index
80!
81! Revision 1.1  1997/08/11 06:26:19  raasch
82! Initial revision
83!
84!
85! Description:
86! ------------
87! Compute the time step under consideration of the FCL and diffusion criterion.
88!------------------------------------------------------------------------------!
89
90    USE arrays_3d
91    USE control_parameters
92    USE cpulog
93    USE grid_variables
94    USE indices
95    USE interfaces
96    USE pegrid
97    USE statistics
98
99    IMPLICIT NONE
100
101    INTEGER ::  i, j, k
102
103    REAL ::  div, dt_diff, dt_diff_l, dt_u, dt_v, dt_w, percent_change, &
104             u_gtrans_l, value, v_gtrans_l
105
106    REAL, DIMENSION(2)         ::  uv_gtrans, uv_gtrans_l
107    REAL, DIMENSION(nzb+1:nzt) ::  dxyz2_min
108
109
110    CALL cpu_log( log_point(12), 'calculate_timestep', 'start' )
111
112!
113!-- Determine the maxima of the velocity components.
114    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, u, 'abs', &
115                         u_max, u_max_ijk )
116    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, v, 'abs', &
117                         v_max, v_max_ijk )
118    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, w, 'abs', &
119                         w_max, w_max_ijk )
120
121!
122!-- In case maxima of the horizontal velocity components have been found at the
123!-- bottom boundary (k=nzb), the corresponding maximum at level k=1 is chosen
124!-- if the Dirichlet-boundary condition ('mirror') has been set. This is
125!-- necessary, because otherwise in case of Galilei-transform a far too large
126!-- velocity (having the respective opposite sign) would be used for the time
127!-- step determination (almost double the mean flow velocity).
128    IF ( ibc_uv_b == 0 )  THEN
129       IF ( u_max_ijk(1) == nzb )  THEN
130          u_max        = -u_max
131          u_max_ijk(1) = nzb + 1
132       ENDIF
133       IF ( v_max_ijk(1) == nzb )  THEN
134          v_max        = -v_max
135          v_max_ijk(1) = nzb + 1
136       ENDIF
137    ENDIF
138
139!
140!-- In case of Galilei-transform not using the geostrophic wind as translation
141!-- velocity, compute the volume-averaged horizontal velocity components, which
142!-- will then be subtracted from the horizontal wind for the time step and
143!-- horizontal advection routines.
144    IF ( galilei_transformation  .AND. .NOT.  use_ug_for_galilei_tr )  THEN
145       IF ( flow_statistics_called )  THEN
146!
147!--       Horizontal averages already existent, just need to average them
148!--       vertically.
149          u_gtrans = 0.0
150          v_gtrans = 0.0
151          DO  k = nzb+1, nzt
152             u_gtrans = u_gtrans + hom(k,1,1,0)
153             v_gtrans = v_gtrans + hom(k,1,2,0)
154          ENDDO
155          u_gtrans = u_gtrans / REAL( nzt - nzb )
156          v_gtrans = v_gtrans / REAL( nzt - nzb )
157       ELSE
158!
159!--       Averaging over the entire model domain.
160          uv_gtrans_l = 0.0
161          DO  i = nxl, nxr
162             DO  j = nys, nyn
163                DO  k = nzb+1, nzt
164                   uv_gtrans_l(1) = uv_gtrans_l(1) + u(k,j,i)
165                   uv_gtrans_l(2) = uv_gtrans_l(2) + v(k,j,i)
166                ENDDO
167             ENDDO
168          ENDDO
169          uv_gtrans_l = uv_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) )
170#if defined( __parallel )
171          CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM, &
172                              comm2d, ierr )
173          u_gtrans = uv_gtrans(1) / REAL( numprocs )
174          v_gtrans = uv_gtrans(2) / REAL( numprocs )
175#else
176          u_gtrans = uv_gtrans_l(1)
177          v_gtrans = uv_gtrans_l(2)
178#endif
179       ENDIF
180    ENDIF
181
182    IF ( .NOT. dt_fixed )  THEN
183!
184!--    Variable time step:
185!
186!--    For each component, compute the maximum time step according to the
187!--    FCL-criterion.
188       dt_u = dx / ( ABS( u_max - u_gtrans ) + 1.0E-10 )
189       dt_v = dy / ( ABS( v_max - v_gtrans ) + 1.0E-10 )
190       dt_w = dzu(MAX( 1, w_max_ijk(1) )) / ( ABS( w_max ) + 1.0E-10 )
191
192!
193!--    Compute time step according to the diffusion criterion.
194!--    First calculate minimum grid spacing which only depends on index k
195!--    Note: also at k=nzb+1 a full grid length is being assumed, although
196!--          in the Prandtl-layer friction term only dz/2 is used.
197!--          Experience from the old model seems to justify this.
198       dt_diff_l = 999999.0
199
200       DO  k = nzb+1, nzt
201           dxyz2_min(k) = MIN( dx2, dy2, dzu(k)*dzu(k) ) * 0.125
202       ENDDO
203
204!$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l)
205!$OMP DO
206       DO  i = nxl, nxr
207          DO  j = nys, nyn
208             DO  k = nzb+1, nzt
209                value = dxyz2_min(k) / ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 )
210
211                dt_diff_l = MIN( value, dt_diff_l )
212             ENDDO
213          ENDDO
214       ENDDO
215!$OMP END PARALLEL
216#if defined( __parallel )
217       CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, &
218                           ierr )
219#else
220       dt_diff = dt_diff_l
221#endif
222
223!
224!--    In case of non-cyclic lateral boundaries, the diffusion time step
225!--    may be further restricted by the lateral damping layer (damping only
226!--    along x and y)
227       IF ( bc_lr /= 'cyclic' )  THEN
228          dt_diff = MIN( dt_diff, 0.125 * dx2 / ( km_damp_max + 1E-20 ) )
229       ELSEIF ( bc_ns /= 'cyclic' )  THEN
230          dt_diff = MIN( dt_diff, 0.125 * dy2 / ( km_damp_max + 1E-20 ) )
231       ENDIF
232
233!
234!--    The time step is the minimum of the 3 components and the diffusion time
235!--    step minus a reduction to be on the safe side. Factor 0.5 is necessary
236!--    since the leap-frog scheme always progresses by 2 * delta t.
237!--    The user has to set the cfl_factor small enough to ensure that the
238!--    divergences do not become too large.
239!--    The time step must not exceed the maximum allowed value.
240       IF ( timestep_scheme(1:5) == 'runge' )  THEN
241          dt_3d = cfl_factor * MIN( dt_diff, dt_u, dt_v, dt_w )
242       ELSE
243          dt_3d = cfl_factor * 0.5 * MIN( dt_diff, dt_u, dt_v, dt_w )
244       ENDIF
245       dt_3d = MIN( dt_3d, dt_max )
246
247!
248!--    Remember the restricting time step criterion for later output.
249       IF ( dt_diff > MIN( dt_u, dt_v, dt_w ) )  THEN
250          timestep_reason = 'A'
251       ELSE
252          timestep_reason = 'D'
253       ENDIF
254
255!
256!--    Set flag if the time step becomes too small.
257       IF ( dt_3d < ( 0.00001 * dt_max ) )  THEN
258          stop_dt = .TRUE.
259          IF ( myid == 0 )  THEN
260             PRINT*,'+++ time_step: Time step has reached minimum limit.'
261             PRINT*,'    dt      = ', dt_3d, ' s   Simulation is terminated.'
262             PRINT*,'    old_dt  = ', old_dt, ' s'
263             PRINT*,'    dt_u    = ', dt_u, ' s'
264             PRINT*,'    dt_v    = ', dt_v, ' s'
265             PRINT*,'    dt_w    = ', dt_w, ' s'
266             PRINT*,'    dt_diff = ', dt_diff, ' s'
267             PRINT*,'    u_max   = ', u_max, ' m/s   k=', u_max_ijk(1), &
268                    '  j=', u_max_ijk(2), '  i=', u_max_ijk(3)
269             PRINT*,'    v_max   = ', v_max, ' m/s   k=', v_max_ijk(1), &
270                    '  j=', v_max_ijk(2), '  i=', v_max_ijk(3)
271             PRINT*,'    w_max   = ', w_max, ' m/s   k=', w_max_ijk(1), &
272                    '  j=', w_max_ijk(2), '  i=', w_max_ijk(3)
273          ENDIF
274       ENDIF
275
276!
277!--    Ensure a smooth value (two significant digits) of the timestep. For
278!--    other schemes than Runge-Kutta, the following restrictions appear:
279!--    The current timestep is only then changed, if the change relative to
280!--    its previous value exceeds +5 % or -2 %. In case of a timestep
281!--    reduction, at least 30 iterations have to be performed before a timestep
282!--    enlargement is permitted again.
283       percent_change = dt_3d / old_dt - 1.0
284       IF ( percent_change > 0.05  .OR.  percent_change < -0.02  .OR. &
285            timestep_scheme(1:5) == 'runge' )  THEN
286
287!
288!--       Time step enlargement by no more than 2 %.
289          IF ( percent_change > 0.0  .AND.  simulated_time /= 0.0  .AND. &
290               timestep_scheme(1:5) /= 'runge' )  THEN
291             dt_3d = 1.02 * old_dt
292          ENDIF
293
294!
295!--       A relatively smooth value of the time step is ensured by taking
296!--       only the first two significant digits.
297          div = 1000.0
298          DO  WHILE ( dt_3d < div )
299             div = div / 10.0
300          ENDDO
301          dt_3d = NINT( dt_3d * 100.0 / div ) * div / 100.0
302
303!
304!--       Now the time step can be adjusted.
305          IF ( percent_change < 0.0  .OR.  timestep_scheme(1:5) == 'runge' ) &
306          THEN
307!
308!--          Time step reduction.
309             old_dt = dt_3d
310             dt_changed = .TRUE.
311          ELSE
312!
313!--          For other timestep schemes , the time step is only enlarged
314!--          after at least 30 iterations since the previous time step
315!--          change or, of course, after model initialization.
316             IF ( current_timestep_number >= last_dt_change + 30  .OR. &
317                  simulated_time == 0.0 )  THEN
318                old_dt = dt_3d
319                dt_changed = .TRUE.
320             ELSE
321                dt_3d = old_dt
322                dt_changed = .FALSE.
323             ENDIF
324
325          ENDIF
326       ELSE
327!
328!--       No time step change since the difference is too small.
329          dt_3d = old_dt
330          dt_changed = .FALSE.
331       ENDIF
332
333       IF ( dt_changed )  last_dt_change = current_timestep_number
334
335    ENDIF
336
337    CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' )
338
339
340 END SUBROUTINE timestep
Note: See TracBrowser for help on using the repository browser.