source: palm/trunk/SOURCE/subsidence_mod.f90 @ 4262

Last change on this file since 4262 was 4182, checked in by scharf, 5 years ago
  • corrected "Former revisions" section
  • minor formatting in "Former revisions" section
  • added "Author" section
  • Property svn:keywords set to Id
File size: 13.8 KB
RevLine 
[1850]1!> @file subsidence_mod.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!
[411]20! Current revisions:
21! -----------------
[1683]22!
[2233]23!
[1383]24! Former revisions:
25! -----------------
26! $Id: subsidence_mod.f90 4182 2019-08-22 15:20:23Z schwenkel $
[4182]27! Corrected "Former revisions" section
28!
29! 3655 2019-01-07 16:51:22Z knoop
[3554]30! add subroutine and variable description
[1383]31!
[4182]32! Revision 3.7 2009-12-11 14:15:58Z heinze
33! Initial revision
34!
[411]35! Description:
36! ------------
[1682]37!> Impact of large-scale subsidence or ascent as tendency term for use
38!> in the prognostic equation of potential temperature. This enables the
39!> construction of a constant boundary layer height z_i with time.
[411]40!-----------------------------------------------------------------------------!
[1682]41 MODULE subsidence_mod
42 
[411]43
44
45    IMPLICIT NONE
46
47    PRIVATE
48    PUBLIC  init_w_subsidence, subsidence
49
50    INTERFACE init_w_subsidence
51       MODULE PROCEDURE init_w_subsidence
52    END INTERFACE init_w_subsidence
53
54    INTERFACE subsidence
55       MODULE PROCEDURE subsidence
56       MODULE PROCEDURE subsidence_ij
57    END INTERFACE subsidence
58
59 CONTAINS
60
[1682]61!------------------------------------------------------------------------------!
62! Description:
63! ------------
[3554]64!> Initialize vertical subsidence velocity w_subs.
[1682]65!------------------------------------------------------------------------------!
[411]66    SUBROUTINE init_w_subsidence 
67
[1320]68       USE arrays_3d,                                                          &
69           ONLY:  dzu, w_subs, zu
[411]70
[1320]71       USE control_parameters,                                                 &
[3294]72           ONLY:  message_string, ocean_mode, subs_vertical_gradient,          &
[1320]73                  subs_vertical_gradient_level, subs_vertical_gradient_level_i
74
75       USE indices,                                                            &
76           ONLY:  nzb, nzt
77
78       USE kinds
79
[411]80       IMPLICIT NONE
81
[3554]82       INTEGER(iwp) ::  i !< loop index
83       INTEGER(iwp) ::  k !< loop index
[411]84
[3554]85       REAL(wp)     ::  gradient   !< vertical gradient of subsidence velocity
86       REAL(wp)     ::  ws_surface !< subsidence velocity at the surface
[1320]87
[1862]88       IF ( .NOT. ALLOCATED( w_subs ) )  THEN
[411]89          ALLOCATE( w_subs(nzb:nzt+1) )
[1353]90          w_subs = 0.0_wp
[411]91       ENDIF
92
[3294]93       IF ( ocean_mode )  THEN
[3302]94          message_string = 'applying large scale vertical motion is not ' //   &
[3294]95                           'allowed for ocean mode'
[411]96          CALL message( 'init_w_subsidence', 'PA0324', 2, 2, 0, 6, 0 )
97       ENDIF
98
99!
100!--   Compute the profile of the subsidence/ascent velocity
101!--   using the given gradients
102      i = 1
[1353]103      gradient = 0.0_wp
104      ws_surface = 0.0_wp
[411]105     
106
[580]107      subs_vertical_gradient_level_i(1) = 0
[411]108      DO  k = 1, nzt+1
[1365]109         IF ( i < 11 )  THEN
[580]110            IF ( subs_vertical_gradient_level(i) < zu(k)  .AND. &
[1353]111                 subs_vertical_gradient_level(i) >= 0.0_wp )  THEN
112               gradient = subs_vertical_gradient(i) / 100.0_wp
[580]113               subs_vertical_gradient_level_i(i) = k - 1
[411]114               i = i + 1
115            ENDIF
116         ENDIF
[1353]117         IF ( gradient /= 0.0_wp )  THEN
[411]118            IF ( k /= 1 )  THEN
119               w_subs(k) = w_subs(k-1) + dzu(k) * gradient
120            ELSE
[1353]121               w_subs(k) = ws_surface   + 0.5_wp * dzu(k) * gradient
[411]122            ENDIF
123         ELSE
124            w_subs(k) = w_subs(k-1)
125         ENDIF
126      ENDDO
127
128!
129!--   In case of no given gradients for the subsidence/ascent velocity,
130!--   choose zero gradient
[1353]131      IF ( subs_vertical_gradient_level(1) == -9999999.9_wp )  THEN
132         subs_vertical_gradient_level(1) = 0.0_wp
[411]133      ENDIF
134
135    END SUBROUTINE init_w_subsidence
136
137
[1682]138!------------------------------------------------------------------------------!
139! Description:
140! ------------
[3554]141!> Add effect of large-scale subsidence to variable.
[1682]142!------------------------------------------------------------------------------!
[1365]143    SUBROUTINE subsidence( tendency, var, var_init, ls_index ) 
[411]144
[1320]145       USE arrays_3d,                                                          &
146           ONLY:  ddzu, w_subs
[411]147
[1320]148       USE control_parameters,                                                 &
[1380]149           ONLY:  dt_3d, intermediate_timestep_count, large_scale_forcing,     &
150                  scalar_rayleigh_damping
[1320]151
152       USE indices,                                                            &
[2232]153           ONLY:  nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt,        &
154                  wall_flags_0
[1320]155
156       USE kinds
157
[1365]158       USE statistics,                                                         &
[1382]159           ONLY:  sums_ls_l, weight_substep
[1365]160
[411]161       IMPLICIT NONE
162 
[3554]163       INTEGER(iwp) ::  i        !< loop index
164       INTEGER(iwp) ::  j        !< loop index
165       INTEGER(iwp) ::  k        !< loop index
166       INTEGER(iwp) ::  ls_index !< index of large-scale subsidence in sums_ls_l
[411]167
[3554]168       REAL(wp)     ::  tmp_tend !< temporary tendency
169       REAL(wp)     ::  tmp_grad !< temporary gradient
[411]170   
[3554]171       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  var      !< variable where to add subsidence
172       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  tendency !< tendency of var
173       REAL(wp), DIMENSION(nzb:nzt+1) ::  var_init                     !< initialization profile of var
174       REAL(wp), DIMENSION(nzb:nzt+1) ::  var_mod                      !< modified profile of var
[411]175
[464]176       var_mod = var_init
[411]177
178!
179!--    Influence of w_subsidence on the current tendency term
180       DO  i = nxl, nxr
181          DO  j = nys, nyn
[1382]182
[2232]183             DO  k = nzb+1, nzt 
[1365]184                IF ( w_subs(k) < 0.0_wp )  THEN    ! large-scale subsidence
185                   tmp_tend = - w_subs(k) *                                    &
[2232]186                              ( var(k+1,j,i) - var(k,j,i) ) * ddzu(k+1) *      &
187                                        MERGE( 1.0_wp, 0.0_wp,                 &
188                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]189                ELSE                               ! large-scale ascent
190                   tmp_tend = - w_subs(k) *                                    &
[2232]191                              ( var(k,j,i) - var(k-1,j,i) ) * ddzu(k) *        &
192                                        MERGE( 1.0_wp, 0.0_wp,                 &
193                                               BTEST( wall_flags_0(k,j,i), 0 ) )
[411]194                ENDIF
[1365]195
196                tendency(k,j,i) = tendency(k,j,i) + tmp_tend
197
198                IF ( large_scale_forcing )  THEN
199                   sums_ls_l(k,ls_index) = sums_ls_l(k,ls_index) + tmp_tend    &
[2232]200                                 * weight_substep(intermediate_timestep_count) &
201                                 * MERGE( 1.0_wp, 0.0_wp,                      &
202                                          BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]203                ENDIF
[411]204             ENDDO
[1382]205
[1862]206             IF ( large_scale_forcing )  THEN
207                sums_ls_l(nzt+1,ls_index) = sums_ls_l(nzt,ls_index)
208             ENDIF
[1382]209
[411]210          ENDDO
211       ENDDO
212
213!
214!--    Shifting of the initial profile is especially necessary with Rayleigh
215!--    damping switched on
[1729]216       IF ( scalar_rayleigh_damping .AND.                                      &
217            intermediate_timestep_count == 1 )  THEN
[1380]218          DO  k = nzb, nzt
219             IF ( w_subs(k) < 0.0_wp )  THEN      ! large-scale subsidence
220                var_mod(k) = var_init(k) - dt_3d * w_subs(k) *  &
221                                  ( var_init(k+1) - var_init(k) ) * ddzu(k+1)
222             ENDIF
223          ENDDO
[411]224!
[1380]225!--      At the upper boundary, the initial profile is shifted with aid of
226!--      the gradient tmp_grad. (This is ok if the gradients are linear.)
227         IF ( w_subs(nzt) < 0.0_wp )  THEN
228            tmp_grad = ( var_init(nzt+1) - var_init(nzt) ) * ddzu(nzt+1)
229            var_mod(nzt+1) = var_init(nzt+1) -  &
230                                 dt_3d * w_subs(nzt+1) * tmp_grad
231         ENDIF
[411]232       
233
[1380]234         DO  k = nzt+1, nzb+1, -1
235            IF ( w_subs(k) >= 0.0_wp )  THEN  ! large-scale ascent
236               var_mod(k) = var_init(k) - dt_3d * w_subs(k) *  &
237                                  ( var_init(k) - var_init(k-1) ) * ddzu(k) 
238            ENDIF
239         ENDDO
240!
241!--      At the lower boundary shifting is not necessary because the
242!--      subsidence velocity w_subs(nzb) vanishes.
243         IF ( w_subs(nzb+1) >= 0.0_wp )  THEN
244            var_mod(nzb) = var_init(nzb)
[411]245         ENDIF
[1380]246
247         var_init = var_mod
[411]248      ENDIF
249
250
251 END SUBROUTINE subsidence
252
[1682]253!------------------------------------------------------------------------------!
254! Description:
255! ------------
[3554]256!> Add effect of large-scale subsidence to variable.
[1682]257!------------------------------------------------------------------------------!
[1365]258 SUBROUTINE subsidence_ij( i, j, tendency, var, var_init, ls_index ) 
[411]259
[1320]260       USE arrays_3d,                                                          &
261           ONLY:  ddzu, w_subs
[411]262
[1320]263       USE control_parameters,                                                 &
[1380]264           ONLY:  dt_3d, intermediate_timestep_count, large_scale_forcing,     &
265                  scalar_rayleigh_damping
[1320]266
267       USE indices,                                                            &
[2232]268           ONLY:  nxl, nxlg, nxrg, nyng, nys, nysg, nzb, nzt, wall_flags_0
[1320]269
270       USE kinds
271
[1365]272       USE statistics,                                                         &
[1382]273           ONLY:  sums_ls_l, weight_substep
[1365]274
[411]275       IMPLICIT NONE
276 
[3554]277       INTEGER(iwp) ::  i        !< loop variable
278       INTEGER(iwp) ::  j        !< loop variable
279       INTEGER(iwp) ::  k        !< loop variable
280       INTEGER(iwp) ::  ls_index !< index of large-scale subsidence in sums_ls_l
[411]281
[3554]282       REAL(wp)     ::  tmp_tend !< temporary tendency
283       REAL(wp)     ::  tmp_grad !< temporary gradient
[411]284   
[3554]285       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  var      !< variable where to add subsidence
286       REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  tendency !< tendency of var
287       REAL(wp), DIMENSION(nzb:nzt+1) ::  var_init                     !< initialization profile of var
288       REAL(wp), DIMENSION(nzb:nzt+1) ::  var_mod                      !< modified profile of var
[411]289
[464]290       var_mod = var_init
[411]291
292!
293!--    Influence of w_subsidence on the current tendency term
[2232]294       DO  k = nzb+1, nzt 
[1365]295          IF ( w_subs(k) < 0.0_wp )  THEN      ! large-scale subsidence
[2232]296             tmp_tend = - w_subs(k) * ( var(k+1,j,i) - var(k,j,i) )            &
297                                    * ddzu(k+1)                                &
298                                    * MERGE( 1.0_wp, 0.0_wp,                   &
299                                             BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]300          ELSE                                 ! large-scale ascent
[2232]301             tmp_tend = - w_subs(k) * ( var(k,j,i) - var(k-1,j,i) ) * ddzu(k)  &
302                                    * MERGE( 1.0_wp, 0.0_wp,                   &
303                                             BTEST( wall_flags_0(k,j,i), 0 ) )
[411]304          ENDIF
[1365]305
306          tendency(k,j,i) = tendency(k,j,i) + tmp_tend
307
308          IF ( large_scale_forcing )  THEN
309             sums_ls_l(k,ls_index) = sums_ls_l(k,ls_index) + tmp_tend          &
[2232]310                                  * weight_substep(intermediate_timestep_count)&
311                                  * MERGE( 1.0_wp, 0.0_wp,                     &
312                                           BTEST( wall_flags_0(k,j,i), 0 ) )
[1365]313          ENDIF
[411]314       ENDDO
315
[1489]316       IF ( large_scale_forcing )  THEN
317          sums_ls_l(nzt+1,ls_index) = sums_ls_l(nzt,ls_index)
318       ENDIF
[411]319
320!
321!--    Shifting of the initial profile is especially necessary with Rayleigh
322!--    damping switched on
[1729]323       IF ( scalar_rayleigh_damping .AND.                                      &
324            intermediate_timestep_count == 1 )  THEN
[1380]325          IF ( i == nxl .AND. j == nys )  THEN ! shifting only once per PE
[411]326
[1380]327             DO  k = nzb, nzt
328                IF ( w_subs(k) < 0.0_wp )  THEN      ! large-scale subsidence
329                   var_mod(k) = var_init(k) - dt_3d * w_subs(k) *  &
330                                     ( var_init(k+1) - var_init(k) ) * ddzu(k+1)
331                ENDIF
332             ENDDO
333!
334!--          At the upper boundary, the initial profile is shifted with aid of
335!--          the gradient tmp_grad. (This is ok if the gradients are linear.)
336             IF ( w_subs(nzt) < 0.0_wp )  THEN
337                tmp_grad = ( var_init(nzt+1) - var_init(nzt) ) * ddzu(nzt+1)
338                var_mod(nzt+1) = var_init(nzt+1) -  &
339                                     dt_3d * w_subs(nzt+1) * tmp_grad
[411]340             ENDIF
341       
342
[1380]343             DO  k = nzt+1, nzb+1, -1
344                IF ( w_subs(k) >= 0.0_wp )  THEN  ! large-scale ascent
345                   var_mod(k) = var_init(k) - dt_3d * w_subs(k) *  &
346                                      ( var_init(k) - var_init(k-1) ) * ddzu(k)
347                ENDIF
348             ENDDO
349!
350!--          At the lower boundary shifting is not necessary because the
351!--          subsidence velocity w_subs(nzb) vanishes.
352             IF ( w_subs(nzb+1) >= 0.0_wp )  THEN
353                var_mod(nzb) = var_init(nzb)
[411]354             ENDIF
355
[1380]356             var_init = var_mod 
[411]357
[1380]358          ENDIF
[411]359       ENDIF
360
361 END SUBROUTINE subsidence_ij
362
363
364 END MODULE subsidence_mod
Note: See TracBrowser for help on using the repository browser.