[1873] | 1 | !> @file advec_s_pw.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 | ! |
---|
[484] | 20 | ! Current revisions: |
---|
[1] | 21 | ! ----------------- |
---|
[1354] | 22 | ! |
---|
[2233] | 23 | ! |
---|
[1321] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: advec_s_pw.f90 3927 2019-04-23 13:24:29Z knoop $ |
---|
[3927] | 27 | ! pointer attribute removed from scalar 3d-array for performance reasons |
---|
| 28 | ! |
---|
| 29 | ! 3665 2019-01-10 08:28:24Z raasch |
---|
[3665] | 30 | ! unused variables removed |
---|
| 31 | ! |
---|
| 32 | ! 3655 2019-01-07 16:51:22Z knoop |
---|
[3636] | 33 | ! nopointer option removed |
---|
| 34 | ! |
---|
| 35 | ! 3547 2018-11-21 13:21:24Z suehring |
---|
[3547] | 36 | ! variables documented |
---|
| 37 | ! |
---|
| 38 | ! 3538 2018-11-20 10:55:41Z suehring |
---|
[3538] | 39 | ! Remove unnecessary double-masking of topography |
---|
| 40 | ! |
---|
| 41 | ! 3302 2018-10-03 02:39:40Z raasch |
---|
[3302] | 42 | ! Stokes drift velocity added |
---|
| 43 | ! |
---|
| 44 | ! 2718 2018-01-02 08:49:38Z maronga |
---|
[2716] | 45 | ! Corrected "Former revisions" section |
---|
| 46 | ! |
---|
| 47 | ! 2696 2017-12-14 17:12:51Z kanani |
---|
| 48 | ! Change in file header (GPL part) |
---|
[1321] | 49 | ! |
---|
[2716] | 50 | ! 2233 2017-05-30 18:08:54Z suehring |
---|
| 51 | ! |
---|
[2233] | 52 | ! 2232 2017-05-30 17:47:52Z suehring |
---|
| 53 | ! topography representation via flags |
---|
| 54 | ! |
---|
[2001] | 55 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
| 56 | ! Forced header and separation lines into 80 columns |
---|
| 57 | ! |
---|
[1874] | 58 | ! 1873 2016-04-18 14:50:06Z maronga |
---|
| 59 | ! Module renamed (removed _mod) |
---|
| 60 | ! |
---|
| 61 | ! |
---|
[1851] | 62 | ! 1850 2016-04-08 13:29:27Z maronga |
---|
| 63 | ! Module renamed |
---|
| 64 | ! |
---|
| 65 | ! |
---|
[1683] | 66 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 67 | ! Code annotations made doxygen readable |
---|
| 68 | ! |
---|
[1375] | 69 | ! 1374 2014-04-25 12:55:07Z raasch |
---|
| 70 | ! missing variables added to ONLY list |
---|
| 71 | ! |
---|
[1354] | 72 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
| 73 | ! REAL constants provided with KIND-attribute |
---|
| 74 | ! |
---|
[1321] | 75 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 76 | ! ONLY-attribute added to USE-statements, |
---|
| 77 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 78 | ! kinds are defined in new module kinds, |
---|
| 79 | ! revision history before 2012 removed, |
---|
| 80 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 81 | ! all variable declaration statements |
---|
[1] | 82 | ! |
---|
[1037] | 83 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 84 | ! code put under GPL (PALM 3.9) |
---|
| 85 | ! |
---|
[1011] | 86 | ! 1010 2012-09-20 07:59:54Z raasch |
---|
| 87 | ! cpp switch __nopointer added for pointer free version |
---|
| 88 | ! |
---|
[1] | 89 | ! Revision 1.1 1997/08/29 08:54:20 raasch |
---|
| 90 | ! Initial revision |
---|
| 91 | ! |
---|
| 92 | ! |
---|
| 93 | ! Description: |
---|
| 94 | ! ------------ |
---|
[1682] | 95 | !> Advection term for scalar variables using the Piacsek and Williams scheme |
---|
| 96 | !> (form C3). Contrary to PW itself, for reasons of accuracy their scheme is |
---|
| 97 | !> slightly modified as follows: the values of those scalars that are used for |
---|
| 98 | !> the computation of the flux divergence are reduced by the value of the |
---|
| 99 | !> relevant scalar at the location where the difference is computed (sk(k,j,i)). |
---|
| 100 | !> NOTE: at the first grid point above the surface computation still takes place! |
---|
[1] | 101 | !------------------------------------------------------------------------------! |
---|
[1682] | 102 | MODULE advec_s_pw_mod |
---|
| 103 | |
---|
[1] | 104 | |
---|
| 105 | PRIVATE |
---|
| 106 | PUBLIC advec_s_pw |
---|
| 107 | |
---|
| 108 | INTERFACE advec_s_pw |
---|
| 109 | MODULE PROCEDURE advec_s_pw |
---|
| 110 | MODULE PROCEDURE advec_s_pw_ij |
---|
| 111 | END INTERFACE |
---|
| 112 | |
---|
| 113 | CONTAINS |
---|
| 114 | |
---|
| 115 | |
---|
| 116 | !------------------------------------------------------------------------------! |
---|
[1682] | 117 | ! Description: |
---|
| 118 | ! ------------ |
---|
| 119 | !> Call for all grid points |
---|
[1] | 120 | !------------------------------------------------------------------------------! |
---|
| 121 | SUBROUTINE advec_s_pw( sk ) |
---|
| 122 | |
---|
[1320] | 123 | USE arrays_3d, & |
---|
[3302] | 124 | ONLY: dd2zu, tend, u, u_stokes_zu, v, v_stokes_zu, w |
---|
[1] | 125 | |
---|
[1320] | 126 | USE control_parameters, & |
---|
| 127 | ONLY: u_gtrans, v_gtrans |
---|
| 128 | |
---|
| 129 | USE grid_variables, & |
---|
| 130 | ONLY: ddx, ddy |
---|
| 131 | |
---|
| 132 | USE indices, & |
---|
[3927] | 133 | ONLY: nxl, nxlg, nxr, nxrg, nyn, nyng, nys, nysg, nzb, nzt |
---|
[1320] | 134 | |
---|
| 135 | USE kinds |
---|
| 136 | |
---|
| 137 | |
---|
[1] | 138 | IMPLICIT NONE |
---|
| 139 | |
---|
[3547] | 140 | INTEGER(iwp) :: i !< grid index along x-direction |
---|
| 141 | INTEGER(iwp) :: j !< grid index along y-direction |
---|
| 142 | INTEGER(iwp) :: k !< grid index along z-direction |
---|
[1] | 143 | |
---|
[3302] | 144 | REAL(wp) :: gu !< local additional advective velocity |
---|
| 145 | REAL(wp) :: gv !< local additional advective velocity |
---|
| 146 | |
---|
[3927] | 147 | REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: sk |
---|
[1] | 148 | |
---|
| 149 | |
---|
| 150 | DO i = nxl, nxr |
---|
| 151 | DO j = nys, nyn |
---|
[2232] | 152 | DO k = nzb+1, nzt |
---|
[3302] | 153 | |
---|
| 154 | ! |
---|
| 155 | !-- Galilean transformation + Stokes drift velocity (ocean only) |
---|
| 156 | gu = u_gtrans - u_stokes_zu(k) |
---|
| 157 | gv = v_gtrans - v_stokes_zu(k) |
---|
| 158 | |
---|
| 159 | tend(k,j,i) = tend(k,j,i) + & |
---|
| 160 | ( -0.5_wp * ( ( u(k,j,i+1) - gu ) * & |
---|
| 161 | ( sk(k,j,i+1) - sk(k,j,i) ) & |
---|
| 162 | - ( u(k,j,i) - gu ) * & |
---|
| 163 | ( sk(k,j,i-1) - sk(k,j,i) ) & |
---|
| 164 | ) * ddx & |
---|
| 165 | -0.5_wp * ( ( v(k,j+1,i) - gv ) * & |
---|
| 166 | ( sk(k,j+1,i) - sk(k,j,i) ) & |
---|
| 167 | - ( v(k,j,i) - gv ) * & |
---|
| 168 | ( sk(k,j-1,i) - sk(k,j,i) ) & |
---|
| 169 | ) * ddy & |
---|
| 170 | - ( w(k,j,i) * & |
---|
| 171 | ( sk(k+1,j,i) - sk(k,j,i) ) & |
---|
| 172 | - w(k-1,j,i) * & |
---|
| 173 | ( sk(k-1,j,i) - sk(k,j,i) ) & |
---|
| 174 | ) * dd2zu(k) & |
---|
[3538] | 175 | ) |
---|
[1] | 176 | ENDDO |
---|
| 177 | ENDDO |
---|
| 178 | ENDDO |
---|
| 179 | |
---|
| 180 | END SUBROUTINE advec_s_pw |
---|
| 181 | |
---|
| 182 | |
---|
| 183 | !------------------------------------------------------------------------------! |
---|
[1682] | 184 | ! Description: |
---|
| 185 | ! ------------ |
---|
| 186 | !> Call for grid point i,j |
---|
[1] | 187 | !------------------------------------------------------------------------------! |
---|
| 188 | SUBROUTINE advec_s_pw_ij( i, j, sk ) |
---|
| 189 | |
---|
[1320] | 190 | USE arrays_3d, & |
---|
[3302] | 191 | ONLY: dd2zu, tend, u, u_stokes_zu, v, v_stokes_zu, w |
---|
[1] | 192 | |
---|
[1320] | 193 | USE control_parameters, & |
---|
| 194 | ONLY: u_gtrans, v_gtrans |
---|
| 195 | |
---|
| 196 | USE grid_variables, & |
---|
| 197 | ONLY: ddx, ddy |
---|
| 198 | |
---|
| 199 | USE indices, & |
---|
[3927] | 200 | ONLY: nxlg, nxrg, nyng, nysg, nzb, nzt |
---|
[1320] | 201 | |
---|
| 202 | USE kinds |
---|
| 203 | |
---|
| 204 | |
---|
[1] | 205 | IMPLICIT NONE |
---|
| 206 | |
---|
[3547] | 207 | INTEGER(iwp) :: i !< grid index along x-direction |
---|
| 208 | INTEGER(iwp) :: j !< grid index along y-direction |
---|
| 209 | INTEGER(iwp) :: k !< grid index along z-direction |
---|
[1] | 210 | |
---|
[3302] | 211 | REAL(wp) :: gu !< local additional advective velocity |
---|
| 212 | REAL(wp) :: gv !< local additional advective velocity |
---|
| 213 | |
---|
[3927] | 214 | REAL(wp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) :: sk |
---|
[1] | 215 | |
---|
| 216 | |
---|
[2232] | 217 | DO k = nzb+1, nzt |
---|
[3302] | 218 | |
---|
| 219 | ! |
---|
| 220 | !-- Galilean transformation + Stokes drift velocity (ocean only) |
---|
| 221 | gu = u_gtrans - u_stokes_zu(k) |
---|
| 222 | gv = v_gtrans - v_stokes_zu(k) |
---|
| 223 | |
---|
| 224 | tend(k,j,i) = tend(k,j,i) + & |
---|
| 225 | ( -0.5_wp * ( ( u(k,j,i+1) - gu ) * & |
---|
| 226 | ( sk(k,j,i+1) - sk(k,j,i) ) & |
---|
| 227 | - ( u(k,j,i) - gu ) * & |
---|
| 228 | ( sk(k,j,i-1) - sk(k,j,i) ) & |
---|
| 229 | ) * ddx & |
---|
| 230 | -0.5_wp * ( ( v(k,j+1,i) - gv ) * & |
---|
| 231 | ( sk(k,j+1,i) - sk(k,j,i) ) & |
---|
| 232 | - ( v(k,j,i) - gv ) * & |
---|
| 233 | ( sk(k,j-1,i) - sk(k,j,i) ) & |
---|
| 234 | ) * ddy & |
---|
| 235 | - ( w(k,j,i) * & |
---|
| 236 | ( sk(k+1,j,i) - sk(k,j,i) ) & |
---|
| 237 | - w(k-1,j,i) * & |
---|
| 238 | ( sk(k-1,j,i) - sk(k,j,i) ) & |
---|
| 239 | ) * dd2zu(k) & |
---|
[3538] | 240 | ) |
---|
[1] | 241 | ENDDO |
---|
| 242 | |
---|
| 243 | END SUBROUTINE advec_s_pw_ij |
---|
| 244 | |
---|
| 245 | END MODULE advec_s_pw_mod |
---|