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

Last change on this file since 4493 was 4444, checked in by raasch, 5 years ago

bugfix: cpp-directives for serial mode added

  • Property svn:keywords set to Id
File size: 15.5 KB
Line 
1!> @file timestep.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: timestep.f90 4444 2020-03-05 15:59:50Z pavelkrc $
27! bugfix: cpp-directives for serial mode added
28!
29! 4360 2020-01-07 11:25:50Z suehring
30! Added missing OpenMP directives
31!
32! 4233 2019-09-20 09:55:54Z knoop
33! OpenACC data update host removed
34!
35! 4182 2019-08-22 15:20:23Z scharf
36! Corrected "Former revisions" section
37!
38! 4101 2019-07-17 15:14:26Z gronemeier
39! - consider 2*Km within diffusion criterion as Km is considered twice within
40!   the diffusion of e,
41! - in RANS mode, instead of considering each wind component individually use
42!   the wind speed of 3d wind vector in CFL criterion
43! - do not limit the increase of dt based on its previous value in RANS mode
44!
45! 3658 2019-01-07 20:28:54Z knoop
46! OpenACC port for SPEC
47!
48! Revision 1.1  1997/08/11 06:26:19  raasch
49! Initial revision
50!
51!
52! Description:
53! ------------
54!> Compute the time step under consideration of the FCL and diffusion criterion.
55!------------------------------------------------------------------------------!
56 SUBROUTINE timestep
57 
58
59    USE arrays_3d,                                                             &
60        ONLY:  dzu, dzw, kh, km, u, u_stokes_zu, v, v_stokes_zu, w
61
62    USE control_parameters,                                                    &
63        ONLY:  cfl_factor, dt_3d, dt_fixed, dt_max, galilei_transformation,    &
64               message_string, rans_mode, stop_dt, timestep_reason, u_gtrans,  &
65               use_ug_for_galilei_tr, v_gtrans
66
67#if defined( __parallel )
68    USE control_parameters,                                                    &
69        ONLY:  coupling_mode, terminate_coupled, terminate_coupled_remote
70#endif
71
72    USE cpulog,                                                                &
73        ONLY:  cpu_log, log_point
74
75    USE grid_variables,                                                        &
76        ONLY:  dx, dx2, dy, dy2
77
78    USE indices,                                                               &
79        ONLY:  nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt
80
81    USE interfaces
82
83    USE kinds
84
85    USE bulk_cloud_model_mod,                                                  &
86        ONLY:  dt_precipitation
87
88    USE pegrid
89
90    USE pmc_interface,                                                         &
91        ONLY:  nested_run
92
93    USE statistics,                                                            &
94        ONLY:  flow_statistics_called, hom, u_max, u_max_ijk, v_max, v_max_ijk,&
95               w_max, w_max_ijk
96
97#if defined( __parallel )
98    USE vertical_nesting_mod,                                                  &
99        ONLY:  vnested, vnest_timestep_sync
100#endif
101
102    IMPLICIT NONE
103
104    INTEGER(iwp) ::  i !<
105    INTEGER(iwp) ::  j !<
106    INTEGER(iwp) ::  k !<
107    INTEGER(iwp) ::  km_max_ijk(3) = -1  !< index values (i,j,k) of location where km_max occurs
108    INTEGER(iwp) ::  kh_max_ijk(3) = -1  !< index values (i,j,k) of location where kh_max occurs
109
110    LOGICAL ::  stop_dt_local !< local switch for controlling the time stepping
111
112    REAL(wp) ::  div               !<
113    REAL(wp) ::  dt_diff           !<
114    REAL(wp) ::  dt_diff_l         !<
115    REAL(wp) ::  dt_u              !<
116    REAL(wp) ::  dt_u_l            !<
117    REAL(wp) ::  dt_v              !<
118    REAL(wp) ::  dt_v_l            !<
119    REAL(wp) ::  dt_w              !<
120    REAL(wp) ::  dt_w_l            !<
121    REAL(wp) ::  km_max            !< maximum of Km in entire domain
122    REAL(wp) ::  kh_max            !< maximum of Kh in entire domain
123    REAL(wp) ::  u_gtrans_l        !<
124    REAL(wp) ::  v_gtrans_l        !<
125 
126    REAL(wp), DIMENSION(2)         ::  uv_gtrans_l !<
127#if defined( __parallel )
128    REAL(wp), DIMENSION(2)         ::  uv_gtrans   !<
129    REAL(wp), DIMENSION(3)         ::  reduce      !<
130    REAL(wp), DIMENSION(3)         ::  reduce_l    !<
131#endif
132    REAL(wp), DIMENSION(nzb+1:nzt) ::  dxyz2_min   !< 
133    !$ACC DECLARE CREATE(dxyz2_min)
134
135
136    CALL cpu_log( log_point(12), 'calculate_timestep', 'start' )
137
138!
139!-- In case of Galilei-transform not using the geostrophic wind as translation
140!-- velocity, compute the volume-averaged horizontal velocity components, which
141!-- will then be subtracted from the horizontal wind for the time step and
142!-- horizontal advection routines.
143    IF ( galilei_transformation  .AND. .NOT.  use_ug_for_galilei_tr )  THEN
144       IF ( flow_statistics_called )  THEN
145!
146!--       Horizontal averages already existent, just need to average them
147!--       vertically.
148          u_gtrans = 0.0_wp
149          v_gtrans = 0.0_wp
150          DO  k = nzb+1, nzt
151             u_gtrans = u_gtrans + hom(k,1,1,0)
152             v_gtrans = v_gtrans + hom(k,1,2,0)
153          ENDDO
154          u_gtrans = u_gtrans / REAL( nzt - nzb, KIND=wp )
155          v_gtrans = v_gtrans / REAL( nzt - nzb, KIND=wp )
156       ELSE
157!
158!--       Averaging over the entire model domain.
159          u_gtrans_l = 0.0_wp
160          v_gtrans_l = 0.0_wp
161          DO  i = nxl, nxr
162             DO  j = nys, nyn
163                DO  k = nzb+1, nzt
164                   u_gtrans_l = u_gtrans_l + u(k,j,i)
165                   v_gtrans_l = v_gtrans_l + v(k,j,i)
166                ENDDO
167             ENDDO
168          ENDDO
169          uv_gtrans_l(1) = u_gtrans_l /                                        &
170                           REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb), KIND=wp )
171          uv_gtrans_l(2) = v_gtrans_l /                                        &
172                           REAL( (nxr-nxl+1)*(nyn-nys+1)*(nzt-nzb), KIND=wp )
173#if defined( __parallel )
174          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
175          CALL MPI_ALLREDUCE( uv_gtrans_l, uv_gtrans, 2, MPI_REAL, MPI_SUM,    &
176                              comm2d, ierr )
177          u_gtrans = uv_gtrans(1) / REAL( numprocs, KIND=wp )
178          v_gtrans = uv_gtrans(2) / REAL( numprocs, KIND=wp )
179#else
180          u_gtrans = uv_gtrans_l(1)
181          v_gtrans = uv_gtrans_l(2)
182#endif
183       ENDIF
184    ENDIF
185
186!
187!-- Determine the maxima of the velocity components, including their
188!-- grid index positions.
189    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, u, 'abs', 0.0_wp, &
190                         u_max, u_max_ijk )
191    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, v, 'abs', 0.0_wp, &
192                         v_max, v_max_ijk )
193    CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, w, 'abs', 0.0_wp, &
194                         w_max, w_max_ijk )
195
196    IF ( .NOT. dt_fixed )  THEN
197!
198!--    Variable time step:
199!--    Calculate the maximum time step according to the CFL-criterion
200       dt_u_l = 999999.9_wp
201       dt_v_l = 999999.9_wp
202       dt_w_l = 999999.9_wp
203
204       IF ( .NOT. rans_mode )  THEN
205!
206!--       Consider each velocity component individually
207
208          !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i,j,k) &
209          !$ACC COPY(dt_u_l, dt_v_l, dt_w_l, u_stokes_zu, v_stokes_zu) &
210          !$ACC REDUCTION(MIN: dt_u_l, dt_v_l, dt_w_l) &
211          !$ACC PRESENT(u, v, w, dzu)
212          !$OMP PARALLEL DO PRIVATE(i,j,k) &
213          !$OMP REDUCTION(MIN: dt_u_l, dt_v_l, dt_w_l)
214          DO  i = nxl, nxr
215             DO  j = nys, nyn
216                DO  k = nzb+1, nzt
217                   dt_u_l = MIN( dt_u_l, ( dx     /                               &
218                                    ( ABS( u(k,j,i) - u_gtrans + u_stokes_zu(k) ) &
219                                      + 1.0E-10_wp ) ) )
220                   dt_v_l = MIN( dt_v_l, ( dy     /                               &
221                                    ( ABS( v(k,j,i) - v_gtrans + v_stokes_zu(k) ) &
222                                      + 1.0E-10_wp ) ) )
223                   dt_w_l = MIN( dt_w_l, ( dzu(k) /                               &
224                                    ( ABS( w(k,j,i) )            + 1.0E-10_wp ) ) )
225                ENDDO
226             ENDDO
227          ENDDO
228
229       ELSE
230!
231!--       Consider the wind speed at the scalar-grid point
232!--       !> @note considering the wind speed instead of each individual wind
233!--       !>       component is only a workaround so far. This might has to be
234!--       !>       changed in the future.
235
236          !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i,j,k) &
237          !$ACC COPY(dt_u_l, u_stokes_zu, v_stokes_zu) &
238          !$ACC REDUCTION(MIN: dt_u_l) &
239          !$ACC PRESENT(u, v, w, dzu)
240          !$OMP PARALLEL DO PRIVATE(i,j,k) &
241          !$OMP REDUCTION(MIN: dt_u_l)
242          DO  i = nxl, nxr
243             DO  j = nys, nyn
244                DO  k = nzb+1, nzt
245                   dt_u_l = MIN( dt_u_l, ( MIN( dx, dy, dzu(k) ) / ( &
246                      SQRT(  ( 0.5 * ( u(k,j,i) + u(k,j,i+1) ) - u_gtrans + u_stokes_zu(k) )**2   &
247                           + ( 0.5 * ( v(k,j,i) + v(k,j+1,i) ) - v_gtrans + v_stokes_zu(k) )**2   &
248                           + ( 0.5 * ( w(k,j,i) + w(k-1,j,i) )                             )**2 ) &
249                      + 1.0E-10_wp ) ) )
250                ENDDO
251             ENDDO
252          ENDDO
253         
254          dt_v_l = dt_u_l
255          dt_w_l = dt_u_l
256
257       ENDIF
258
259#if defined( __parallel )
260       reduce_l(1) = dt_u_l
261       reduce_l(2) = dt_v_l
262       reduce_l(3) = dt_w_l
263       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
264       CALL MPI_ALLREDUCE( reduce_l, reduce, 3, MPI_REAL, MPI_MIN, comm2d, ierr )
265       dt_u = reduce(1)
266       dt_v = reduce(2)
267       dt_w = reduce(3)
268#else
269       dt_u = dt_u_l
270       dt_v = dt_v_l
271       dt_w = dt_w_l
272#endif
273
274!
275!--    Compute time step according to the diffusion criterion.
276!--    First calculate minimum grid spacing which only depends on index k.
277!--    When using the dynamic subgrid model, negative km are possible.
278       dt_diff_l = 999999.0_wp
279
280       !$ACC PARALLEL LOOP PRESENT(dxyz2_min, dzw)
281       DO  k = nzb+1, nzt
282           dxyz2_min(k) = MIN( dx2, dy2, dzw(k)*dzw(k) ) * 0.125_wp
283       ENDDO
284
285       !$OMP PARALLEL private(i,j,k) reduction(MIN: dt_diff_l)
286       !$OMP DO
287       !$ACC PARALLEL LOOP COLLAPSE(3) PRIVATE(i,j,k) &
288       !$ACC COPY(dt_diff_l) REDUCTION(MIN: dt_diff_l) &
289       !$ACC PRESENT(dxyz2_min, kh, km)
290       DO  i = nxl, nxr
291          DO  j = nys, nyn
292             DO  k = nzb+1, nzt
293                dt_diff_l = MIN( dt_diff_l,                                       &
294                                 dxyz2_min(k) /                                   &
295                                    ( MAX( kh(k,j,i), 2.0_wp * ABS( km(k,j,i) ) ) &
296                                      + 1E-20_wp ) )
297             ENDDO
298          ENDDO
299       ENDDO
300       !$OMP END PARALLEL
301#if defined( __parallel )
302       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
303       CALL MPI_ALLREDUCE( dt_diff_l, dt_diff, 1, MPI_REAL, MPI_MIN, comm2d,   &
304                           ierr )
305#else
306       dt_diff = dt_diff_l
307#endif
308
309!
310!--    The time step is the minimum of the 3-4 components and the diffusion time
311!--    step minus a reduction (cfl_factor) to be on the safe side.
312!--    The time step must not exceed the maximum allowed value.
313       dt_3d = cfl_factor * MIN( dt_diff, dt_u, dt_v, dt_w, dt_precipitation )
314       dt_3d = MIN( dt_3d, dt_max )
315
316!
317!--    Remember the restricting time step criterion for later output.
318       IF ( MIN( dt_u, dt_v, dt_w ) < dt_diff )  THEN
319          timestep_reason = 'A'
320       ELSE
321          timestep_reason = 'D'
322       ENDIF
323
324!
325!--    Set flag if the time step becomes too small.
326       IF ( dt_3d < ( 0.00001_wp * dt_max ) )  THEN
327          stop_dt = .TRUE.
328
329!
330!--       Determine the maxima of the diffusion coefficients, including their
331!--       grid index positions.
332          CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, km, 'abs',  &
333                               0.0_wp, km_max, km_max_ijk )
334          CALL global_min_max( nzb, nzt+1, nysg, nyng, nxlg, nxrg, kh, 'abs',  &
335                               0.0_wp, kh_max, kh_max_ijk )
336
337          WRITE( message_string, * ) 'Time step has reached minimum limit.',   &
338               '&dt              = ', dt_3d, ' s  Simulation is terminated.',  &
339               '&dt_u            = ', dt_u, ' s',                              &
340               '&dt_v            = ', dt_v, ' s',                              &
341               '&dt_w            = ', dt_w, ' s',                              &
342               '&dt_diff         = ', dt_diff, ' s',                           &
343               '&u_max           = ', u_max, ' m/s    k=', u_max_ijk(1),       &
344               '  j=', u_max_ijk(2), '  i=', u_max_ijk(3),                     &
345               '&v_max           = ', v_max, ' m/s    k=', v_max_ijk(1),       &
346               '  j=', v_max_ijk(2), '  i=', v_max_ijk(3),                     &
347               '&w_max           = ', w_max, ' m/s    k=', w_max_ijk(1),       &
348               '  j=', w_max_ijk(2), '  i=', w_max_ijk(3),                     &
349               '&km_max          = ', km_max, ' m2/s2  k=', km_max_ijk(1),     &
350               '  j=', km_max_ijk(2), '  i=', km_max_ijk(3),                   &
351               '&kh_max          = ', kh_max, ' m2/s2  k=', kh_max_ijk(1),     &
352                '  j=', kh_max_ijk(2), '  i=', kh_max_ijk(3)
353          CALL message( 'timestep', 'PA0312', 0, 1, 0, 6, 0 )
354!
355!--       In case of coupled runs inform the remote model of the termination
356!--       and its reason, provided the remote model has not already been
357!--       informed of another termination reason (terminate_coupled > 0) before.
358#if defined( __parallel )
359          IF ( coupling_mode /= 'uncoupled' .AND. terminate_coupled == 0 )  THEN
360             terminate_coupled = 2
361             IF ( myid == 0 )  THEN
362                CALL MPI_SENDRECV( &
363                     terminate_coupled,        1, MPI_INTEGER, target_id,  0,  &
364                     terminate_coupled_remote, 1, MPI_INTEGER, target_id,  0,  &
365                     comm_inter, status, ierr )
366             ENDIF
367             CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0,      &
368                             comm2d, ierr)
369          ENDIF
370#endif
371       ENDIF
372
373!
374!--    In case of nested runs all parent/child processes have to terminate if
375!--    one process has set the stop flag, i.e. they need to set the stop flag
376!--    too.
377       IF ( nested_run )  THEN
378          stop_dt_local = stop_dt
379#if defined( __parallel )
380          CALL MPI_ALLREDUCE( stop_dt_local, stop_dt, 1, MPI_LOGICAL, MPI_LOR, &
381                              MPI_COMM_WORLD, ierr )
382#endif
383       ENDIF
384
385!
386!--    Ensure a smooth value (two significant digits) of the timestep.
387       div = 1000.0_wp
388       DO  WHILE ( dt_3d < div )
389          div = div / 10.0_wp
390       ENDDO
391       dt_3d = NINT( dt_3d * 100.0_wp / div ) * div / 100.0_wp
392
393    ENDIF
394
395#if defined( __parallel )
396!
397!-- Vertical nesting: coarse and fine grid timestep has to be identical   
398    IF ( vnested )  CALL vnest_timestep_sync
399#endif
400
401    CALL cpu_log( log_point(12), 'calculate_timestep', 'stop' )
402
403 END SUBROUTINE timestep
Note: See TracBrowser for help on using the repository browser.