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

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