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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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