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

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