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

Last change on this file since 255 was 226, checked in by raasch, 16 years ago

preparations for the next release

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