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

Last change on this file since 203 was 110, checked in by raasch, 17 years ago

New:
---
Allows runs for a coupled atmosphere-ocean LES,
coupling frequency is controlled by new d3par-parameter dt_coupling,
the coupling mode (atmosphere_to_ocean or ocean_to_atmosphere) for the
respective processes is read from environment variable coupling_mode,
which is set by the mpiexec-command,
communication between the two models is done using the intercommunicator
comm_inter,
local files opened by the ocean model get the additional suffic "_O".
Assume saturation at k=nzb_s_inner(j,i) for atmosphere coupled to ocean.

A momentum flux can be set as top boundary condition using the new
inipar parameter top_momentumflux_u|v.

Non-cyclic boundary conditions can be used along all horizontal directions.

Quantities w*p* and w"e can be output as vertical profiles.

Initial profiles are reset to constant profiles in case that initializing_actions /= 'set_constant_profiles'. (init_rankine)

Optionally calculate km and kh from initial TKE e_init.

Changed:


Remaining variables iran changed to iran_part (advec_particles, init_particles).

In case that the presure solver is not called for every Runge-Kutta substep
(call_psolver_at_all_substeps = .F.), it is called after the first substep
instead of the last. In that case, random perturbations are also added to the
velocity field after the first substep.

Initialization of km,kh = 0.00001 for ocean = .T. (for ocean = .F. it remains 0.01).

Allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of cloud_physics = .T.).

Errors:


Bugs from code parts for non-cyclic boundary conditions are removed: loops for
u and v are starting from index nxlu, nysv, respectively. The radiation boundary
condition is used for every Runge-Kutta substep. Velocity phase speeds for
the radiation boundary conditions are calculated for the first Runge-Kutta
substep only and reused for the further substeps. New arrays c_u, c_v, and c_w
are defined for this purpose. Several index errors are removed from the
radiation boundary condition code parts. Upper bounds for calculating
u_0 and v_0 (in production_e) are nxr+1 and nyn+1 because otherwise these
values are not available in case of non-cyclic boundary conditions.

+dots_num_palm in module user, +module netcdf_control in user_init (both in user_interface)

Bugfix: wrong sign removed from the buoyancy production term in the case use_reference = .T. (production_e)

Bugfix: Error message concerning output of particle concentration (pc) modified (check_parameters).

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