source: palm/trunk/SOURCE/diffusion_w.f90 @ 3436

Last change on this file since 3436 was 3241, checked in by raasch, 6 years ago

various changes to avoid compiler warnings (mainly removal of unused variables)

  • Property svn:keywords set to Id
File size: 19.8 KB
RevLine 
[1873]1!> @file diffusion_w.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!
[2718]17! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[484]20! Current revisions:
[1]21! -----------------
[1341]22!
[2233]23!
[1321]24! Former revisions:
25! -----------------
26! $Id: diffusion_w.f90 3241 2018-09-12 15:02:00Z gronemeier $
[3241]27! unused variables removed
28!
29! 2718 2018-01-02 08:49:38Z maronga
[2716]30! Corrected "Former revisions" section
31!
32! 2696 2017-12-14 17:12:51Z kanani
33! Change in file header (GPL part)
[1321]34!
[2716]35! 2233 2017-05-30 18:08:54Z suehring
36!
[2233]37! 2232 2017-05-30 17:47:52Z suehring
38! Adjustments to new topography and surface concept
39!
[2119]40! 2118 2017-01-17 16:38:49Z raasch
41! OpenACC version of subroutine removed
42!
[2038]43! 2037 2016-10-26 11:15:40Z knoop
44! Anelastic approximation implemented
45!
[2001]46! 2000 2016-08-20 18:09:15Z knoop
47! Forced header and separation lines into 80 columns
48!
[1874]49! 1873 2016-04-18 14:50:06Z maronga
50! Module renamed (removed _mod)
51!
[1851]52! 1850 2016-04-08 13:29:27Z maronga
53! Module renamed
54!
[1683]55! 1682 2015-10-07 23:56:08Z knoop
56! Code annotations made doxygen readable
57!
[1375]58! 1374 2014-04-25 12:55:07Z raasch
59! vsws + vswst removed from acc-present-list
60!
[1354]61! 1353 2014-04-08 15:21:23Z heinze
62! REAL constants provided with KIND-attribute
63!
[1341]64! 1340 2014-03-25 19:45:13Z kanani
65! REAL constants defined as wp-kind
66!
[1323]67! 1322 2014-03-20 16:38:49Z raasch
68! REAL constants defined as wp-kind
69!
[1321]70! 1320 2014-03-20 08:40:49Z raasch
[1320]71! ONLY-attribute added to USE-statements,
72! kind-parameters added to all INTEGER and REAL declaration statements,
73! kinds are defined in new module kinds,
74! revision history before 2012 removed,
75! comment fields (!:) to be used for variable explanations added to
76! all variable declaration statements
[1]77!
[1258]78! 1257 2013-11-08 15:18:40Z raasch
79! openacc loop and loop vector clauses removed, declare create moved after
80! the FORTRAN declaration statement
81!
[1132]82! 1128 2013-04-12 06:19:32Z raasch
83! loop index bounds in accelerator version replaced by i_left, i_right, j_south,
84! j_north
85!
[1037]86! 1036 2012-10-22 13:43:42Z raasch
87! code put under GPL (PALM 3.9)
88!
[1017]89! 1015 2012-09-27 09:23:24Z raasch
90! accelerator version (*_acc) added
91!
[1002]92! 1001 2012-09-13 14:08:46Z raasch
93! arrays comunicated by module instead of parameter list
94!
[979]95! 978 2012-08-09 08:28:32Z fricke
96! outflow damping layer removed
97! kmxm_x/_z and kmxp_x/_z change to kmxm and kmxp
98! kmym_y/_z and kmyp_y/_z change to kmym and kmyp
99!
[1]100! Revision 1.1  1997/09/12 06:24:11  raasch
101! Initial revision
102!
103!
104! Description:
105! ------------
[1682]106!> Diffusion term of the w-component
[1]107!------------------------------------------------------------------------------!
[1682]108 MODULE diffusion_w_mod
109 
[1]110
111    PRIVATE
[2118]112    PUBLIC diffusion_w
[1]113
114    INTERFACE diffusion_w
115       MODULE PROCEDURE diffusion_w
116       MODULE PROCEDURE diffusion_w_ij
117    END INTERFACE diffusion_w
118
119 CONTAINS
120
121
122!------------------------------------------------------------------------------!
[1682]123! Description:
124! ------------
125!> Call for all grid points
[1]126!------------------------------------------------------------------------------!
[1001]127    SUBROUTINE diffusion_w
[1]128
[1320]129       USE arrays_3d,                                                          &         
[2037]130           ONLY :  ddzu, ddzw, km, tend, u, v, w, drho_air_zw, rho_air
[1320]131           
132       USE grid_variables,                                                     &     
[2232]133           ONLY :  ddx, ddy
[1320]134           
135       USE indices,                                                            &           
[2232]136           ONLY :  nxl, nxr, nyn, nys, nzb, nzt, wall_flags_0
[1320]137           
138       USE kinds
[1]139
[2232]140       USE surface_mod,                                                        &
141           ONLY :  surf_def_v, surf_lsm_v, surf_usm_v
142
[1]143       IMPLICIT NONE
144
[2232]145       INTEGER(iwp) ::  i             !< running index x direction
146       INTEGER(iwp) ::  j             !< running index y direction
147       INTEGER(iwp) ::  k             !< running index z direction
148       INTEGER(iwp) ::  l             !< running index of surface type, south- or north-facing wall
149       INTEGER(iwp) ::  m             !< running index surface elements
150       INTEGER(iwp) ::  surf_e        !< End index of surface elements at (j,i)-gridpoint
151       INTEGER(iwp) ::  surf_s        !< Start index of surface elements at (j,i)-gridpoint
[1320]152       
[2232]153       REAL(wp) ::  flag              !< flag to mask topography grid points
154       REAL(wp) ::  kmxm              !<
155       REAL(wp) ::  kmxp              !<
156       REAL(wp) ::  kmym              !<
157       REAL(wp) ::  kmyp              !<
158       REAL(wp) ::  mask_west         !< flag to mask vertical wall west of the grid point
159       REAL(wp) ::  mask_east         !< flag to mask vertical wall east of the grid point
160       REAL(wp) ::  mask_south        !< flag to mask vertical wall south of the grid point
161       REAL(wp) ::  mask_north        !< flag to mask vertical wall north of the grid point
[1001]162
[1]163
164
165       DO  i = nxl, nxr
166          DO  j = nys, nyn
[2232]167             DO  k = nzb+1, nzt-1
[1]168!
[2232]169!--             Predetermine flag to mask topography and wall-bounded grid points.
170                flag       = MERGE( 1.0_wp, 0.0_wp,                            &
171                                    BTEST( wall_flags_0(k,j,i),   3 ) ) 
172                mask_east  = MERGE( 1.0_wp, 0.0_wp,                            &
173                                    BTEST( wall_flags_0(k,j,i+1), 3 ) )
174                mask_west  = MERGE( 1.0_wp, 0.0_wp,                            &
175                                    BTEST( wall_flags_0(k,j,i-1), 3 ) )
176                mask_south = MERGE( 1.0_wp, 0.0_wp,                            &
177                                    BTEST( wall_flags_0(k,j-1,i), 3 ) )
178                mask_north = MERGE( 1.0_wp, 0.0_wp,                            &
179                                    BTEST( wall_flags_0(k,j+1,i), 3 ) )
180!
[1]181!--             Interpolate eddy diffusivities on staggered gridpoints
[2232]182                kmxp = 0.25_wp * ( km(k,j,i)   +   km(k,j,i+1) +               &
183                                   km(k+1,j,i) + km(k+1,j,i+1) )
184                kmxm = 0.25_wp * ( km(k,j,i)   + km(k,j,i-1)   +               &
185                                   km(k+1,j,i) + km(k+1,j,i-1) )
186                kmyp = 0.25_wp * ( km(k,j,i)   + km(k+1,j,i)   +               &
187                                   km(k,j+1,i) + km(k+1,j+1,i) )
188                kmym = 0.25_wp * ( km(k,j,i)   + km(k+1,j,i)   +               &
189                                   km(k,j-1,i) + km(k+1,j-1,i) )
[1]190
191                tend(k,j,i) = tend(k,j,i)                                      &
[2232]192                       + ( mask_east *  kmxp * (                               &
193                                   ( w(k,j,i+1)   - w(k,j,i)   ) * ddx         &
194                                 + ( u(k+1,j,i+1) - u(k,j,i+1) ) * ddzu(k+1)   &
195                                               )                               &
196                         - mask_west * kmxm *  (                               &
197                                   ( w(k,j,i)     - w(k,j,i-1) ) * ddx         &
198                                 + ( u(k+1,j,i)   - u(k,j,i)   ) * ddzu(k+1)   &
199                                               )                               &
200                         ) * ddx                                 * flag        &
201                       + ( mask_north * kmyp * (                               &
202                                   ( w(k,j+1,i)   - w(k,j,i)   ) * ddy         &
203                                 + ( v(k+1,j+1,i) - v(k,j+1,i) ) * ddzu(k+1)   &
204                                               )                               &
205                         - mask_south * kmym * (                               &
206                                   ( w(k,j,i)     - w(k,j-1,i) ) * ddy         &
207                                 + ( v(k+1,j,i)   - v(k,j,i)   ) * ddzu(k+1)   &
208                                               )                               &
209                         ) * ddy                                 * flag        &
210                       + 2.0_wp * (                                            &
211                         km(k+1,j,i) * ( w(k+1,j,i) - w(k,j,i) )   * ddzw(k+1) &
212                                     * rho_air(k+1)                            &
213                       - km(k,j,i)   * ( w(k,j,i)   - w(k-1,j,i) ) * ddzw(k)   &
214                                     * rho_air(k)                              &
215                                  ) * ddzu(k+1) * drho_air_zw(k) * flag
[1]216             ENDDO
217
218!
[2232]219!--          Add horizontal momentum flux v'w' at north- (l=0) and south-facing (l=1)
220!--          surfaces. Note, in the the flat case, loops won't be entered as
221!--          start_index > end_index. Furtermore, note, no vertical natural surfaces
222!--          so far.           
223!--          Default-type surfaces
224             DO  l = 0, 1
225                surf_s = surf_def_v(l)%start_index(j,i)
226                surf_e = surf_def_v(l)%end_index(j,i)
227                DO  m = surf_s, surf_e
228                   k           = surf_def_v(l)%k(m)
229                   tend(k,j,i) = tend(k,j,i) +                                 &
230                                     surf_def_v(l)%mom_flux_w(m) * ddy
231                ENDDO   
232             ENDDO
[1]233!
[2232]234!--          Natural-type surfaces
235             DO  l = 0, 1
236                surf_s = surf_lsm_v(l)%start_index(j,i)
237                surf_e = surf_lsm_v(l)%end_index(j,i)
238                DO  m = surf_s, surf_e
239                   k           = surf_lsm_v(l)%k(m)
240                   tend(k,j,i) = tend(k,j,i) +                                 &
241                                     surf_lsm_v(l)%mom_flux_w(m) * ddy
242                ENDDO   
243             ENDDO
244!
245!--          Urban-type surfaces
246             DO  l = 0, 1
247                surf_s = surf_usm_v(l)%start_index(j,i)
248                surf_e = surf_usm_v(l)%end_index(j,i)
249                DO  m = surf_s, surf_e
250                   k           = surf_usm_v(l)%k(m)
251                   tend(k,j,i) = tend(k,j,i) +                                 &
252                                     surf_usm_v(l)%mom_flux_w(m) * ddy
253                ENDDO   
254             ENDDO
255!
256!--          Add horizontal momentum flux u'w' at east- (l=2) and west-facing (l=3)
257!--          surface.
258!--          Default-type surfaces
259             DO  l = 2, 3
260                surf_s = surf_def_v(l)%start_index(j,i)
261                surf_e = surf_def_v(l)%end_index(j,i)
262                DO  m = surf_s, surf_e
263                   k           = surf_def_v(l)%k(m)
264                   tend(k,j,i) = tend(k,j,i) +                                 &
265                                     surf_def_v(l)%mom_flux_w(m) * ddx
266                ENDDO   
267             ENDDO
268!
269!--          Natural-type surfaces
270             DO  l = 2, 3
271                surf_s = surf_lsm_v(l)%start_index(j,i)
272                surf_e = surf_lsm_v(l)%end_index(j,i)
273                DO  m = surf_s, surf_e
274                   k           = surf_lsm_v(l)%k(m)
275                   tend(k,j,i) = tend(k,j,i) +                                 &
276                                     surf_lsm_v(l)%mom_flux_w(m) * ddx
277                ENDDO   
278             ENDDO
279!
280!--          Urban-type surfaces
281             DO  l = 2, 3
282                surf_s = surf_usm_v(l)%start_index(j,i)
283                surf_e = surf_usm_v(l)%end_index(j,i)
284                DO  m = surf_s, surf_e
285                   k           = surf_usm_v(l)%k(m)
286                   tend(k,j,i) = tend(k,j,i) +                                 &
287                                     surf_usm_v(l)%mom_flux_w(m) * ddx
288                ENDDO   
289             ENDDO
[1]290
291          ENDDO
292       ENDDO
293
294    END SUBROUTINE diffusion_w
295
296
297!------------------------------------------------------------------------------!
[1682]298! Description:
299! ------------
300!> Call for grid point i,j
[1]301!------------------------------------------------------------------------------!
[1001]302    SUBROUTINE diffusion_w_ij( i, j )
[1]303
[1320]304       USE arrays_3d,                                                          &         
[2037]305           ONLY :  ddzu, ddzw, km, tend, u, v, w, drho_air_zw, rho_air
[1320]306           
307       USE grid_variables,                                                     &     
[2232]308           ONLY :  ddx, ddy
[1320]309           
310       USE indices,                                                            &           
[3241]311           ONLY :  nzb, nzt, wall_flags_0
[1320]312           
313       USE kinds
[1]314
[2232]315       USE surface_mod,                                                        &
316           ONLY :  surf_def_v, surf_lsm_v, surf_usm_v
317
[1]318       IMPLICIT NONE
319
[2232]320
321       INTEGER(iwp) ::  i             !< running index x direction
322       INTEGER(iwp) ::  j             !< running index y direction
323       INTEGER(iwp) ::  k             !< running index z direction
324       INTEGER(iwp) ::  l             !< running index of surface type, south- or north-facing wall
325       INTEGER(iwp) ::  m             !< running index surface elements
326       INTEGER(iwp) ::  surf_e        !< End index of surface elements at (j,i)-gridpoint
327       INTEGER(iwp) ::  surf_s        !< Start index of surface elements at (j,i)-gridpoint
[1320]328       
[2232]329       REAL(wp) ::  flag              !< flag to mask topography grid points
330       REAL(wp) ::  kmxm              !<
331       REAL(wp) ::  kmxp              !<
332       REAL(wp) ::  kmym              !<
333       REAL(wp) ::  kmyp              !<
334       REAL(wp) ::  mask_west         !< flag to mask vertical wall west of the grid point
335       REAL(wp) ::  mask_east         !< flag to mask vertical wall east of the grid point
336       REAL(wp) ::  mask_south        !< flag to mask vertical wall south of the grid point
337       REAL(wp) ::  mask_north        !< flag to mask vertical wall north of the grid point
[1]338
339
[2232]340       DO  k = nzb+1, nzt-1
[1]341!
[2232]342!--       Predetermine flag to mask topography and wall-bounded grid points.
343          flag       = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i),   3 ) ) 
344          mask_east  = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i+1), 3 ) )
345          mask_west  = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i-1), 3 ) )
346          mask_south = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j-1,i), 3 ) )
347          mask_north = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j+1,i), 3 ) )
348!
[1]349!--       Interpolate eddy diffusivities on staggered gridpoints
[1322]350          kmxp = 0.25_wp * ( km(k,j,i)+km(k,j,i+1)+km(k+1,j,i)+km(k+1,j,i+1) )
351          kmxm = 0.25_wp * ( km(k,j,i)+km(k,j,i-1)+km(k+1,j,i)+km(k+1,j,i-1) )
352          kmyp = 0.25_wp * ( km(k,j,i)+km(k+1,j,i)+km(k,j+1,i)+km(k+1,j+1,i) )
353          kmym = 0.25_wp * ( km(k,j,i)+km(k+1,j,i)+km(k,j-1,i)+km(k+1,j-1,i) )
[1]354
355          tend(k,j,i) = tend(k,j,i)                                            &
[2232]356                       + ( mask_east *  kmxp * (                               &
357                                   ( w(k,j,i+1)   - w(k,j,i)   ) * ddx         &
358                                 + ( u(k+1,j,i+1) - u(k,j,i+1) ) * ddzu(k+1)   &
359                                               )                               &
360                         - mask_west * kmxm *  (                               &
361                                   ( w(k,j,i)     - w(k,j,i-1) ) * ddx         &
362                                 + ( u(k+1,j,i)   - u(k,j,i)   ) * ddzu(k+1)   &
363                                               )                               &
364                         ) * ddx                                 * flag        &
365                       + ( mask_north * kmyp * (                               &
366                                   ( w(k,j+1,i)   - w(k,j,i)   ) * ddy         &
367                                 + ( v(k+1,j+1,i) - v(k,j+1,i) ) * ddzu(k+1)   &
368                                               )                               &
369                         - mask_south * kmym * (                               &
370                                   ( w(k,j,i)     - w(k,j-1,i) ) * ddy         &
371                                 + ( v(k+1,j,i)   - v(k,j,i)   ) * ddzu(k+1)   &
372                                               )                               &
373                         ) * ddy                                 * flag        &
374                       + 2.0_wp * (                                            &
375                         km(k+1,j,i) * ( w(k+1,j,i) - w(k,j,i) ) * ddzw(k+1)   &
376                                     * rho_air(k+1)                            &
377                       - km(k,j,i)   * ( w(k,j,i)   - w(k-1,j,i) ) * ddzw(k)   &
378                                     * rho_air(k)                              &
379                                  ) * ddzu(k+1) * drho_air_zw(k) * flag
[1]380       ENDDO
381!
[2232]382!--    Add horizontal momentum flux v'w' at north- (l=0) and south-facing (l=1)
383!--    surfaces. Note, in the the flat case, loops won't be entered as
384!--    start_index > end_index. Furtermore, note, no vertical natural surfaces
385!--    so far.           
386!--    Default-type surfaces
387       DO  l = 0, 1
388          surf_s = surf_def_v(l)%start_index(j,i)
389          surf_e = surf_def_v(l)%end_index(j,i)
390          DO  m = surf_s, surf_e
391             k           = surf_def_v(l)%k(m)
392             tend(k,j,i) = tend(k,j,i) +                                       &
393                                     surf_def_v(l)%mom_flux_w(m) * ddy
394          ENDDO   
395       ENDDO
[51]396!
[2232]397!--    Natural-type surfaces
398       DO  l = 0, 1
399          surf_s = surf_lsm_v(l)%start_index(j,i)
400          surf_e = surf_lsm_v(l)%end_index(j,i)
401          DO  m = surf_s, surf_e
402             k           = surf_lsm_v(l)%k(m)
403             tend(k,j,i) = tend(k,j,i) +                                       &
404                                     surf_lsm_v(l)%mom_flux_w(m) * ddy
405          ENDDO   
406       ENDDO
[1]407!
[2232]408!--    Urban-type surfaces
409       DO  l = 0, 1
410          surf_s = surf_usm_v(l)%start_index(j,i)
411          surf_e = surf_usm_v(l)%end_index(j,i)
412          DO  m = surf_s, surf_e
413             k           = surf_usm_v(l)%k(m)
414             tend(k,j,i) = tend(k,j,i) +                                       &
415                                     surf_usm_v(l)%mom_flux_w(m) * ddy
416          ENDDO   
417       ENDDO
418!
419!--    Add horizontal momentum flux u'w' at east- (l=2) and west-facing (l=3)
420!--    surfaces.
421!--    Default-type surfaces
422       DO  l = 2, 3
423          surf_s = surf_def_v(l)%start_index(j,i)
424          surf_e = surf_def_v(l)%end_index(j,i)
425          DO  m = surf_s, surf_e
426             k           = surf_def_v(l)%k(m)
427             tend(k,j,i) = tend(k,j,i) +                                       &
428                                     surf_def_v(l)%mom_flux_w(m) * ddx
429          ENDDO   
430       ENDDO
431!
432!--    Natural-type surfaces
433       DO  l = 2, 3
434          surf_s = surf_lsm_v(l)%start_index(j,i)
435          surf_e = surf_lsm_v(l)%end_index(j,i)
436          DO  m = surf_s, surf_e
437             k           = surf_lsm_v(l)%k(m)
438             tend(k,j,i) = tend(k,j,i) +                                       &
439                                     surf_lsm_v(l)%mom_flux_w(m) * ddx
440          ENDDO   
441       ENDDO
442!
443!--    Urban-type surfaces
444       DO  l = 2, 3
445          surf_s = surf_usm_v(l)%start_index(j,i)
446          surf_e = surf_usm_v(l)%end_index(j,i)
447          DO  m = surf_s, surf_e
448             k           = surf_usm_v(l)%k(m)
449             tend(k,j,i) = tend(k,j,i) +                                       &
450                                     surf_usm_v(l)%mom_flux_w(m) * ddx
451          ENDDO   
452       ENDDO
[1]453
454
455    END SUBROUTINE diffusion_w_ij
456
[1323]457 END MODULE diffusion_w_mod
Note: See TracBrowser for help on using the repository browser.