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

Last change on this file since 362 was 343, checked in by maronga, 15 years ago

adjustments for lcxt4 and ibmy, allow user 2d xy cross section output at z=nzb+1

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