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

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

OpenACC support added to global_min_max.f90

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