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

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

last commit documented

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