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

Last change on this file since 275 was 274, checked in by heinze, 16 years ago

Indentation of the message calls corrected

  • Property svn:keywords set to Id
File size: 10.6 KB
Line 
1 SUBROUTINE timestep
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6! Output of messages replaced by message handling routine.
7!
8!
9! Former revisions:
10! -----------------
11! $Id: timestep.f90 274 2009-03-26 15:11:21Z raasch $
12!
13! 222 2009-01-12 16:04:16Z letzel
14! Implementation of a MPI-1 Coupling: replaced myid with target_id
15! Bugfix for nonparallel execution
16!
17! 108 2007-08-24 15:10:38Z letzel
18! modifications to terminate coupled runs
19!
20! RCS Log replace by Id keyword, revision history cleaned up
21!
22! Revision 1.21  2006/02/23 12:59:44  raasch
23! nt_anz renamed current_timestep_number
24!
25! Revision 1.1  1997/08/11 06:26:19  raasch
26! Initial revision
27!
28!
29! Description:
30! ------------
31! Compute the time step under consideration of the FCL and diffusion criterion.
32!------------------------------------------------------------------------------!
33
34    USE arrays_3d
35    USE control_parameters
36    USE cpulog
37    USE grid_variables
38    USE indices
39    USE interfaces
40    USE pegrid
41    USE statistics
42
43    IMPLICIT NONE
44
45    INTEGER ::  i, j, k
46
47    REAL ::  div, dt_diff, dt_diff_l, dt_u, dt_v, dt_w, percent_change, &
48             u_gtrans_l, value, v_gtrans_l
49
50    REAL, DIMENSION(2)         ::  uv_gtrans, uv_gtrans_l
51    REAL, DIMENSION(nzb+1:nzt) ::  dxyz2_min
52
53
54    CALL cpu_log( log_point(12), 'calculate_timestep', 'start' )
55
56!
57!-- Determine the maxima of the velocity components.
58    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, u, 'abs', &
59                         u_max, u_max_ijk )
60    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, v, 'abs', &
61                         v_max, v_max_ijk )
62    CALL global_min_max( nzb, nzt+1, nys-1, nyn+1, nxl-1, nxr+1, w, 'abs', &
63                         w_max, w_max_ijk )
64
65!
66!-- In case maxima of the horizontal velocity components have been found at the
67!-- bottom boundary (k=nzb), the corresponding maximum at level k=1 is chosen
68!-- if the Dirichlet-boundary condition ('mirror') has been set. This is
69!-- necessary, because otherwise in case of Galilei-transform a far too large
70!-- velocity (having the respective opposite sign) would be used for the time
71!-- step determination (almost double the mean flow velocity).
72    IF ( ibc_uv_b == 0 )  THEN
73       IF ( u_max_ijk(1) == nzb )  THEN
74          u_max        = -u_max
75          u_max_ijk(1) = nzb + 1
76       ENDIF
77       IF ( v_max_ijk(1) == nzb )  THEN
78          v_max        = -v_max
79          v_max_ijk(1) = nzb + 1
80       ENDIF
81    ENDIF
82
83!
84!-- In case of Galilei-transform not using the geostrophic wind as translation
85!-- velocity, compute the volume-averaged horizontal velocity components, which
86!-- will then be subtracted from the horizontal wind for the time step and
87!-- horizontal advection routines.
88    IF ( galilei_transformation  .AND. .NOT.  use_ug_for_galilei_tr )  THEN
89       IF ( flow_statistics_called )  THEN
90!
91!--       Horizontal averages already existent, just need to average them
92!--       vertically.
93          u_gtrans = 0.0
94          v_gtrans = 0.0
95          DO  k = nzb+1, nzt
96             u_gtrans = u_gtrans + hom(k,1,1,0)
97             v_gtrans = v_gtrans + hom(k,1,2,0)
98          ENDDO
99          u_gtrans = u_gtrans / REAL( nzt - nzb )
100          v_gtrans = v_gtrans / REAL( nzt - nzb )
101       ELSE
102!
103!--       Averaging over the entire model domain.
104          uv_gtrans_l = 0.0
105          DO  i = nxl, nxr
106             DO  j = nys, nyn
107                DO  k = nzb+1, nzt
108                   uv_gtrans_l(1) = uv_gtrans_l(1) + u(k,j,i)
109                   uv_gtrans_l(2) = uv_gtrans_l(2) + v(k,j,i)
110                ENDDO
111             ENDDO
112          ENDDO
113          uv_gtrans_l = uv_gtrans_l / REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb) )
114#if defined( __parallel )
115          CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM, &
116                              comm2d, ierr )
117          u_gtrans = uv_gtrans(1) / REAL( numprocs )
118          v_gtrans = uv_gtrans(2) / REAL( numprocs )
119#else
120          u_gtrans = uv_gtrans_l(1)
121          v_gtrans = uv_gtrans_l(2)
122#endif
123       ENDIF
124    ENDIF
125
126    IF ( .NOT. dt_fixed )  THEN
127!
128!--    Variable time step:
129!
130!--    For each component, compute the maximum time step according to the
131!--    FCL-criterion.
132       dt_u = dx / ( ABS( u_max - u_gtrans ) + 1.0E-10 )
133       dt_v = dy / ( ABS( v_max - v_gtrans ) + 1.0E-10 )
134       dt_w = dzu(MAX( 1, w_max_ijk(1) )) / ( ABS( w_max ) + 1.0E-10 )
135
136!
137!--    Compute time step according to the diffusion criterion.
138!--    First calculate minimum grid spacing which only depends on index k
139!--    Note: also at k=nzb+1 a full grid length is being assumed, although
140!--          in the Prandtl-layer friction term only dz/2 is used.
141!--          Experience from the old model seems to justify this.
142       dt_diff_l = 999999.0
143
144       DO  k = nzb+1, nzt
145           dxyz2_min(k) = MIN( dx2, dy2, dzu(k)*dzu(k) ) * 0.125
146       ENDDO
147
148!$OMP PARALLEL private(i,j,k,value) reduction(MIN: dt_diff_l)
149!$OMP DO
150       DO  i = nxl, nxr
151          DO  j = nys, nyn
152             DO  k = nzb+1, nzt
153                value = dxyz2_min(k) / ( MAX( kh(k,j,i), km(k,j,i) ) + 1E-20 )
154
155                dt_diff_l = MIN( value, dt_diff_l )
156             ENDDO
157          ENDDO
158       ENDDO
159!$OMP END PARALLEL
160#if defined( __parallel )
161       CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d, &
162                           ierr )
163#else
164       dt_diff = dt_diff_l
165#endif
166
167!
168!--    In case of non-cyclic lateral boundaries, the diffusion time step
169!--    may be further restricted by the lateral damping layer (damping only
170!--    along x and y)
171       IF ( bc_lr /= 'cyclic' )  THEN
172          dt_diff = MIN( dt_diff, 0.125 * dx2 / ( km_damp_max + 1E-20 ) )
173       ELSEIF ( bc_ns /= 'cyclic' )  THEN
174          dt_diff = MIN( dt_diff, 0.125 * dy2 / ( km_damp_max + 1E-20 ) )
175       ENDIF
176
177!
178!--    The time step is the minimum of the 3 components and the diffusion time
179!--    step minus a reduction to be on the safe side. Factor 0.5 is necessary
180!--    since the leap-frog scheme always progresses by 2 * delta t.
181!--    The user has to set the cfl_factor small enough to ensure that the
182!--    divergences do not become too large.
183!--    The time step must not exceed the maximum allowed value.
184       IF ( timestep_scheme(1:5) == 'runge' )  THEN
185          dt_3d = cfl_factor * MIN( dt_diff, dt_u, dt_v, dt_w )
186       ELSE
187          dt_3d = cfl_factor * 0.5 * MIN( dt_diff, dt_u, dt_v, dt_w )
188       ENDIF
189       dt_3d = MIN( dt_3d, dt_max )
190
191!
192!--    Remember the restricting time step criterion for later output.
193       IF ( dt_diff > MIN( dt_u, dt_v, dt_w ) )  THEN
194          timestep_reason = 'A'
195       ELSE
196          timestep_reason = 'D'
197       ENDIF
198
199!
200!--    Set flag if the time step becomes too small.
201       IF ( dt_3d < ( 0.00001 * dt_max ) )  THEN
202          stop_dt = .TRUE.
203
204          WRITE( message_string, * ) 'Time step has reached minimum limit.',   &
205                        '&dt      = ', dt_3d, ' s  Simulation is terminated.', &
206                        '&old_dt  = ', old_dt, ' s',                           &
207                        '&dt_u    = ', dt_u, ' s',                             &
208                        '&dt_v    = ', dt_v, ' s',                             &
209                        '&dt_w    = ', dt_w, ' s',                             &
210                        '&dt_diff = ', dt_diff, ' s',                          &
211                        '&u_max   = ', u_max, ' m/s   k=', u_max_ijk(1),       &
212                        '  j=', u_max_ijk(2), '  i=', u_max_ijk(3),            &
213                        '&v_max   = ', v_max, ' m/s   k=', v_max_ijk(1),       &
214                        '  j=', v_max_ijk(2), '  i=', v_max_ijk(3),            &
215                        '&w_max   = ', w_max, ' m/s   k=', w_max_ijk(1),       &
216                        '  j=', w_max_ijk(2), '  i=', w_max_ijk(3)
217          CALL message( 'timestep', 'PA0312', 0, 1, 0, 6, 0 )
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.