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

Last change on this file since 667 was 667, checked in by suehring, 13 years ago

summary:


Gryschka:

  • Coupling with different resolution and different numbers of PEs in ocean and atmosphere is available
  • Exchange of u and v from ocean surface to atmosphere surface
  • Mirror boundary condition for u and v at the bottom are replaced by dirichlet boundary conditions
  • Inflow turbulence is now defined by flucuations around spanwise mean
  • Bugfixes for cyclic_fill and constant_volume_flow

Suehring:

  • New advection added ( Wicker and Skamarock 5th order ), therefore:
    • New module advec_ws.f90
    • Modified exchange of ghost boundaries.
    • Modified evaluation of turbulent fluxes
    • New index bounds nxlg, nxrg, nysg, nyng

advec_ws.f90


Advection scheme for scalars and momentum using the flux formulation of
Wicker and Skamarock 5th order.
Additionally the module contains of a routine using for initialisation and
steering of the statical evaluation. The computation of turbulent fluxes takes
place inside the advection routines.
In case of vector architectures Dirichlet and Radiation boundary conditions are
outstanding and not available. Furthermore simulations within topography are
not possible so far. A further routine local_diss_ij is available and is used
if a control of dissipative fluxes is desired.

check_parameters.f90


Exchange of parameters between ocean and atmosphere via PE0
Check for illegal combination of ws-scheme and timestep scheme.
Check for topography and ws-scheme.
Check for not cyclic boundary conditions in combination with ws-scheme and
loop_optimization = 'vector'.
Check for call_psolver_at_all_substeps and ws-scheme for momentum_advec.

Different processor/grid topology in atmosphere and ocean is now allowed!
Bugfixes in checking for conserve_volume_flow_mode.

exchange_horiz.f90


Dynamic exchange of ghost points with nbgp_local to ensure that no useless
ghost points exchanged in case of multigrid. type_yz(0) and type_xz(0) used for
normal grid, the remaining types used for the several grid levels.
Exchange is done via MPI-Vectors with a dynamic value of ghost points which
depend on the advection scheme. Exchange of left and right PEs is 10% faster
with MPI-Vectors than without.

flow_statistics.f90


When advection is computed with ws-scheme, turbulent fluxes are already
computed in the respective advection routines and buffered in arrays
sums_xxxx_ws_l(). This is due to a consistent treatment of statistics
with the numerics and to avoid unphysical kinks near the surface. So some if-
requests has to be done to dicern between fluxes from ws-scheme other advection
schemes. Furthermore the computation of z_i is only done if the heat flux
exceeds a minimum value. This affects only simulations of a neutral boundary
layer and is due to reasons of computations in the advection scheme.

inflow_turbulence.f90


Using nbgp recycling planes for a better resolution of the turbulent flow near
the inflow.

init_grid.f90


Definition of new array bounds nxlg, nxrg, nysg, nyng on each PE.
Furthermore the allocation of arrays and steering of loops is done with these
parameters. Call of exchange_horiz are modified.
In case of dirichlet bounday condition at the bottom zu(0)=0.0
dzu_mg has to be set explicitly for a equally spaced grid near bottom.
ddzu_pres added to use a equally spaced grid near bottom.

init_pegrid.f90


Moved determination of target_id's from init_coupling
Determination of parameters needed for coupling (coupling_topology, ngp_a, ngp_o)
with different grid/processor-topology in ocean and atmosphere

Adaption of ngp_xy, ngp_y to a dynamic number of ghost points.
The maximum_grid_level changed from 1 to 0. 0 is the normal grid, 1 to
maximum_grid_level the grids for multigrid, in which 0 and 1 are normal grids.
This distinction is due to reasons of data exchange and performance for the
normal grid and grids in poismg.
The definition of MPI-Vectors adapted to a dynamic numer of ghost points.
New MPI-Vectors for data exchange between left and right boundaries added.
This is due to reasons of performance (10% faster).

ATTENTION: nnz_x undefined problem still has to be solved!!!!!!!!
TEST OUTPUT (TO BE REMOVED) logging mpi2 ierr values

parin.f90


Steering parameter dissipation_control added in inipar.

Makefile


Module advec_ws added.

Modules


Removed u_nzb_p1_for_vfc and v_nzb_p1_for_vfc

For coupling with different resolution in ocean and atmophere:
+nx_a, +nx_o, ny_a, +ny_o, ngp_a, ngp_o, +total_2d_o, +total_2d_a,
+coupling_topology

Buffer arrays for the left sided advective fluxes added in arrays_3d.
+flux_s_u, +flux_s_v, +flux_s_w, +diss_s_u, +diss_s_v, +diss_s_w,
+flux_s_pt, +diss_s_pt, +flux_s_e, +diss_s_e, +flux_s_q, +diss_s_q,
+flux_s_sa, +diss_s_sa
3d arrays for dissipation control added. (only necessary for vector arch.)
+var_x, +var_y, +var_z, +gamma_x, +gamma_y, +gamma_z
Default of momentum_advec and scalar_advec changed to 'ws-scheme' .
+exchange_mg added in control_parameters to steer the data exchange.
Parameters +nbgp, +nxlg, +nxrg, +nysg, +nyng added in indices.
flag array +boundary_flags added in indices to steer the degradation of order
of the advective fluxes when non-cyclic boundaries are used.
MPI-datatypes +type_y, +type_y_int and +type_yz for data_exchange added in
pegrid.
+sums_wsus_ws_l, +sums_wsvs_ws_l, +sums_us2_ws_l, +sums_vs2_ws_l,
+sums_ws2_ws_l, +sums_wspts_ws_l, +sums_wssas_ws_l, +sums_wsqs_ws_l
and +weight_substep added in statistics to steer the statistical evaluation
of turbulent fluxes in the advection routines.
LOGICALS +ws_scheme_sca and +ws_scheme_mom added to get a better performance
in prognostic_equations.
LOGICAL +dissipation_control control added to steer numerical dissipation
in ws-scheme.

Changed length of string run_description_header

pres.f90


New allocation of tend when ws-scheme and multigrid is used. This is due to
reasons of perforance of the data_exchange. The same is done with p after
poismg is called.
nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng when no
multigrid is used. Calls of exchange_horiz are modified.

bugfix: After pressure correction no volume flow correction in case of
non-cyclic boundary conditions
(has to be done only before pressure correction)

Call of SOR routine is referenced with ddzu_pres.

prognostic_equations.f90


Calls of the advection routines with WS5 added.
Calls of ws_statistics added to set the statistical arrays to zero after each
time step.

advec_particles.f90


Declaration of de_dx, de_dy, de_dz adapted to additional ghost points.
Furthermore the calls of exchange_horiz were modified.

asselin_filter.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

average_3d_data.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

boundary_conds.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng
Removed mirror boundary conditions for u and v at the bottom in case of
ibc_uv_b == 0. Instead, dirichelt boundary conditions (u=v=0) are set
in init_3d_model

calc_liquid_water_content.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

calc_spectra.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng for
allocation of tend.

check_open.f90


Output of total array size was adapted to nbgp.

data_output_2d.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and
allocation of arrays local_2d and total_2d.
Calls of exchange_horiz are modified.

data_output_2d.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and
allocation of arrays. Calls of exchange_horiz are modified.
Skip-value skip_do_avs changed to a dynamic adaption of ghost points.

data_output_mask.f90


Calls of exchange_horiz are modified.

diffusion_e.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusion_s.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusion_u.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusion_v.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusion_w.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusivities.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

diffusivities.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.
Calls of exchange_horiz are modified.

exchange_horiz_2d.f90


Dynamic exchange of ghost points with nbgp, which depends on the advection
scheme. Exchange between left and right PEs is now done with MPI-vectors.

global_min_max.f90


Adapting of the index arrays, because MINLOC assumes lowerbound
at 1 and not at nbgp.

init_3d_model.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and
allocation of arrays. Calls of exchange_horiz are modified.
Call ws_init to initialize arrays needed for statistical evaluation and
optimization when ws-scheme is used.
Initial volume flow is now calculated by using the variable hom_sum.
Therefore the correction of initial volume flow for non-flat topography
removed (removed u_nzb_p1_for_vfc and v_nzb_p1_for_vfc)
Changed surface boundary conditions for u and v in case of ibc_uv_b == 0 from
mirror bc to dirichlet boundary conditions (u=v=0), so that k=nzb is
representative for the height z0

Bugfix: type conversion of '1' to 64bit for the MAX function (ngp_3d_inner)

init_coupling.f90


determination of target_id's moved to init_pegrid

init_pt_anomaly.f90


Call of exchange_horiz are modified.

init_rankine.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.
Calls of exchange_horiz are modified.

init_slope.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

header.f90


Output of advection scheme.

poismg.f90


Calls of exchange_horiz are modified.

prandtl_fluxes.f90


Changed surface boundary conditions for u and v from mirror bc to dirichelt bc,
therefore u(uzb,:,:) and v(nzb,:,:) is now representative for the height z0
nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

production_e.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng

read_3d_binary.f90


+/- 1 replaced with +/- nbgp when swapping and allocating variables.

sor.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.
Call of exchange_horiz are modified.
bug removed in declaration of ddzw(), nz replaced by nzt+1

subsidence.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

sum_up_3d_data.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

surface_coupler.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in
MPI_SEND() and MPI_RECV.
additional case for nonequivalent processor and grid topopolgy in ocean and
atmosphere added (coupling_topology = 1)

Added exchange of u and v from Ocean to Atmosphere

time_integration.f90


Calls of exchange_horiz are modified.
Adaption to slooping surface.

timestep.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

user_3d_data_averaging.f90, user_data_output_2d.f90, user_data_output_3d.f90,
user_actions.f90, user_init.f90, user_init_plant_canopy.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

user_read_restart_data.f90


Allocation with nbgp.

wall_fluxes.f90


nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng.

write_compressed.f90


Array bounds and nx, ny adapted with nbgp.

sor.f90


bug removed in declaration of ddzw(), nz replaced by nzt+1

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