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

Last change on this file since 622 was 622, checked in by raasch, 13 years ago

New:
---

Optional barriers included in order to speed up collective operations
MPI_ALLTOALL and MPI_ALLREDUCE. This feature is controlled with new initial
parameter collective_wait. Default is .FALSE, but .TRUE. on SGI-type
systems. (advec_particles, advec_s_bc, buoyancy, check_for_restart,
cpu_statistics, data_output_2d, data_output_ptseries, flow_statistics,
global_min_max, inflow_turbulence, init_3d_model, init_particles, init_pegrid,
init_slope, parin, pres, poismg, set_particle_attributes, timestep,
read_var_list, user_statistics, write_compressed, write_var_list)

Adjustments for Kyushu Univ. (lcrte, ibmku). Concerning hybrid
(MPI/openMP) runs, the number of openMP threads per MPI tasks can now
be given as an argument to mrun-option -O. (mbuild, mrun, subjob)

Changed:


Initialization of the module command changed for SGI-ICE/lcsgi (mbuild, subjob)

Errors:


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