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

Last change on this file since 978 was 978, checked in by fricke, 12 years ago

merge fricke branch back into trunk

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