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

Last change on this file since 4274 was 4237, checked in by knoop, 5 years ago

Added missing OpenMP directives within "disturb_field", "surface_layer_fluxes" and "timestep"

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