[1873] | 1 | !> @file advec_v_pw.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[2696] | 3 | !------------------------------------------------------------------------------! |
---|
| 4 | ! This file is part of the PALM model system. |
---|
[1036] | 5 | ! |
---|
[2000] | 6 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
| 7 | ! terms of the GNU General Public License as published by the Free Software |
---|
| 8 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
| 9 | ! version. |
---|
[1036] | 10 | ! |
---|
| 11 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 12 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 13 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 14 | ! |
---|
| 15 | ! You should have received a copy of the GNU General Public License along with |
---|
| 16 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 17 | ! |
---|
[3655] | 18 | ! Copyright 1997-2019 Leibniz Universitaet Hannover |
---|
[2000] | 19 | !------------------------------------------------------------------------------! |
---|
[1036] | 20 | ! |
---|
[484] | 21 | ! Current revisions: |
---|
[1] | 22 | ! ----------------- |
---|
[1354] | 23 | ! |
---|
[2233] | 24 | ! |
---|
[1321] | 25 | ! Former revisions: |
---|
| 26 | ! ----------------- |
---|
| 27 | ! $Id: advec_v_pw.f90 4182 2019-08-22 15:20:23Z scharf $ |
---|
[4182] | 28 | ! Corrected "Former revisions" section |
---|
| 29 | ! |
---|
| 30 | ! 3655 2019-01-07 16:51:22Z knoop |
---|
[3547] | 31 | ! variables documented |
---|
[1321] | 32 | ! |
---|
[4182] | 33 | ! Revision 1.1 1997/08/11 06:09:57 raasch |
---|
| 34 | ! Initial revision |
---|
| 35 | ! |
---|
| 36 | ! |
---|
[1] | 37 | ! Description: |
---|
| 38 | ! ------------ |
---|
[1682] | 39 | !> Advection term for v velocity-component using Piacsek and Williams. |
---|
| 40 | !> Vertical advection at the first grid point above the surface is done with |
---|
| 41 | !> normal centred differences, because otherwise no information from the surface |
---|
| 42 | !> would be communicated upwards due to w=0 at K=nzb. |
---|
[1] | 43 | !------------------------------------------------------------------------------! |
---|
[1682] | 44 | MODULE advec_v_pw_mod |
---|
| 45 | |
---|
[1] | 46 | |
---|
| 47 | PRIVATE |
---|
| 48 | PUBLIC advec_v_pw |
---|
| 49 | |
---|
| 50 | INTERFACE advec_v_pw |
---|
| 51 | MODULE PROCEDURE advec_v_pw |
---|
| 52 | MODULE PROCEDURE advec_v_pw_ij |
---|
| 53 | END INTERFACE advec_v_pw |
---|
| 54 | |
---|
| 55 | CONTAINS |
---|
| 56 | |
---|
| 57 | |
---|
| 58 | !------------------------------------------------------------------------------! |
---|
[1682] | 59 | ! Description: |
---|
| 60 | ! ------------ |
---|
| 61 | !> Call for all grid points |
---|
[1] | 62 | !------------------------------------------------------------------------------! |
---|
| 63 | SUBROUTINE advec_v_pw |
---|
| 64 | |
---|
[1320] | 65 | USE arrays_3d, & |
---|
| 66 | ONLY: ddzw, tend, u, v, w |
---|
[1] | 67 | |
---|
[1320] | 68 | USE control_parameters, & |
---|
| 69 | ONLY: u_gtrans, v_gtrans |
---|
| 70 | |
---|
| 71 | USE grid_variables, & |
---|
| 72 | ONLY: ddx, ddy |
---|
| 73 | |
---|
| 74 | USE indices, & |
---|
[3538] | 75 | ONLY: nxl, nxr, nyn, nysv, nzb, nzt |
---|
[1320] | 76 | |
---|
| 77 | USE kinds |
---|
| 78 | |
---|
| 79 | |
---|
[1] | 80 | IMPLICIT NONE |
---|
| 81 | |
---|
[3547] | 82 | INTEGER(iwp) :: i !< grid index along x-direction |
---|
| 83 | INTEGER(iwp) :: j !< grid index along y-direction |
---|
| 84 | INTEGER(iwp) :: k !< grid index along z-direction |
---|
[1320] | 85 | |
---|
[3547] | 86 | REAL(wp) :: gu !< Galilei-transformation velocity along x |
---|
| 87 | REAL(wp) :: gv !< Galilei-transformation velocity along y |
---|
[1] | 88 | |
---|
| 89 | |
---|
[1353] | 90 | gu = 2.0_wp * u_gtrans |
---|
| 91 | gv = 2.0_wp * v_gtrans |
---|
[1] | 92 | DO i = nxl, nxr |
---|
[106] | 93 | DO j = nysv, nyn |
---|
[2232] | 94 | DO k = nzb+1, nzt |
---|
[1353] | 95 | tend(k,j,i) = tend(k,j,i) - 0.25_wp * ( & |
---|
[1] | 96 | ( v(k,j,i+1) * ( u(k,j-1,i+1) + u(k,j,i+1) - gu ) & |
---|
| 97 | - v(k,j,i-1) * ( u(k,j-1,i) + u(k,j,i) - gu ) ) * ddx & |
---|
| 98 | + ( v(k,j+1,i) * ( v(k,j+1,i) + v(k,j,i) - gv ) & |
---|
| 99 | - v(k,j-1,i) * ( v(k,j,i) + v(k,j-1,i) - gv ) ) * ddy & |
---|
| 100 | + ( v(k+1,j,i) * ( w(k,j-1,i) + w(k,j,i) ) & |
---|
| 101 | - v(k-1,j,i) * ( w(k-1,j-1,i) + w(k-1,j,i) ) ) & |
---|
| 102 | * ddzw(k) & |
---|
[3538] | 103 | ) |
---|
[1] | 104 | ENDDO |
---|
| 105 | ENDDO |
---|
| 106 | ENDDO |
---|
| 107 | |
---|
| 108 | END SUBROUTINE advec_v_pw |
---|
| 109 | |
---|
| 110 | |
---|
| 111 | !------------------------------------------------------------------------------! |
---|
[1682] | 112 | ! Description: |
---|
| 113 | ! ------------ |
---|
| 114 | !> Call for grid point i,j |
---|
[1] | 115 | !------------------------------------------------------------------------------! |
---|
| 116 | SUBROUTINE advec_v_pw_ij( i, j ) |
---|
| 117 | |
---|
[1320] | 118 | USE arrays_3d, & |
---|
| 119 | ONLY: ddzw, tend, u, v, w |
---|
[1] | 120 | |
---|
[1320] | 121 | USE control_parameters, & |
---|
| 122 | ONLY: u_gtrans, v_gtrans |
---|
| 123 | |
---|
| 124 | USE grid_variables, & |
---|
| 125 | ONLY: ddx, ddy |
---|
| 126 | |
---|
| 127 | USE indices, & |
---|
[3538] | 128 | ONLY: nzb, nzt |
---|
[1320] | 129 | |
---|
| 130 | USE kinds |
---|
| 131 | |
---|
| 132 | |
---|
[1] | 133 | IMPLICIT NONE |
---|
| 134 | |
---|
[3547] | 135 | INTEGER(iwp) :: i !< grid index along x-direction |
---|
| 136 | INTEGER(iwp) :: j !< grid index along y-direction |
---|
| 137 | INTEGER(iwp) :: k !< grid index along z-direction |
---|
[1320] | 138 | |
---|
[3547] | 139 | REAL(wp) :: gu !< Galilei-transformation velocity along x |
---|
| 140 | REAL(wp) :: gv !< Galilei-transformation velocity along y |
---|
[1] | 141 | |
---|
| 142 | |
---|
[1353] | 143 | gu = 2.0_wp * u_gtrans |
---|
| 144 | gv = 2.0_wp * v_gtrans |
---|
[2232] | 145 | DO k = nzb+1, nzt |
---|
[1353] | 146 | tend(k,j,i) = tend(k,j,i) - 0.25_wp * ( & |
---|
[1] | 147 | ( v(k,j,i+1) * ( u(k,j-1,i+1) + u(k,j,i+1) - gu ) & |
---|
| 148 | - v(k,j,i-1) * ( u(k,j-1,i) + u(k,j,i) - gu ) ) * ddx & |
---|
| 149 | + ( v(k,j+1,i) * ( v(k,j+1,i) + v(k,j,i) - gv ) & |
---|
| 150 | - v(k,j-1,i) * ( v(k,j,i) + v(k,j-1,i) - gv ) ) * ddy & |
---|
| 151 | + ( v(k+1,j,i) * ( w(k,j-1,i) + w(k,j,i) ) & |
---|
| 152 | - v(k-1,j,i) * ( w(k-1,j-1,i) + w(k-1,j,i) ) ) & |
---|
| 153 | * ddzw(k) & |
---|
[3538] | 154 | ) |
---|
[1] | 155 | ENDDO |
---|
| 156 | |
---|
| 157 | END SUBROUTINE advec_v_pw_ij |
---|
| 158 | |
---|
| 159 | END MODULE advec_v_pw_mod |
---|
| 160 | |
---|