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

Last change on this file since 108 was 108, checked in by letzel, 17 years ago
  • Improved coupler: evaporation - salinity-flux coupling for humidity = .T.,

avoid MPI hangs when coupled runs terminate, add DOC/app/chapter_3.8;

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