source: palm/trunk/SOURCE/advec_v_up.f90 @ 1849

Last change on this file since 1849 was 1818, checked in by maronga, 9 years ago

last commit documented / copyright update

  • Property svn:keywords set to Id
File size: 7.2 KB
Line 
1!> @file advec_v_up.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2016 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: advec_v_up.f90 1818 2016-04-06 15:53:27Z hoffmann $
26!
27! 1682 2015-10-07 23:56:08Z knoop
28! Code annotations made doxygen readable
29!
30! 1353 2014-04-08 15:21:23Z heinze
31! REAL constants provided with KIND-attribute
32!
33! 1320 2014-03-20 08:40:49Z raasch
34! ONLY-attribute added to USE-statements,
35! kind-parameters added to all INTEGER and REAL declaration statements,
36! kinds are defined in new module kinds,
37! revision history before 2012 removed,
38! comment fields (!:) to be used for variable explanations added to
39! all variable declaration statements
40!
41! 1036 2012-10-22 13:43:42Z raasch
42! code put under GPL (PALM 3.9)
43!
44! Revision 1.1  1997/08/29 08:56:05  raasch
45! Initial revision
46!
47!
48! Description:
49! ------------
50!> Advection term for the v velocity-component using upstream scheme.
51!> NOTE: vertical advection at k=1 still has wrong grid spacing for w>0!
52!>       The same problem occurs for all topography boundaries!
53!------------------------------------------------------------------------------!
54 MODULE advec_v_up_mod
55 
56
57    PRIVATE
58    PUBLIC advec_v_up
59
60    INTERFACE advec_v_up
61       MODULE PROCEDURE advec_v_up
62       MODULE PROCEDURE advec_v_up_ij
63    END INTERFACE advec_v_up
64
65 CONTAINS
66
67
68!------------------------------------------------------------------------------!
69! Description:
70! ------------
71!> Call for all grid points
72!------------------------------------------------------------------------------!
73    SUBROUTINE advec_v_up
74
75       USE arrays_3d,                                                          &
76           ONLY:  ddzu, tend, u, v, w
77
78       USE control_parameters,                                                 &
79           ONLY:  u_gtrans, v_gtrans
80
81       USE grid_variables,                                                     &
82           ONLY:  ddx, ddy
83
84       USE indices,                                                            &
85           ONLY:  nxl, nxr, nyn, nysv, nzb_v_inner, nzt
86
87       USE kinds
88
89
90       IMPLICIT NONE
91
92       INTEGER(iwp) ::  i !<
93       INTEGER(iwp) ::  j !<
94       INTEGER(iwp) ::  k !<
95
96       REAL(wp) ::  ukomp !<
97       REAL(wp) ::  vkomp !<
98       REAL(wp) ::  wkomp !<       
99
100
101       DO  i = nxl, nxr
102          DO  j = nysv, nyn
103             DO  k = nzb_v_inner(j,i)+1, nzt
104!
105!--             x-direction
106                ukomp = 0.25_wp * ( u(k,j,i)   + u(k,j-1,i) +                  &
107                                 u(k,j,i+1) + u(k,j-1,i+1) ) - u_gtrans
108                IF ( ukomp > 0.0_wp )  THEN
109                   tend(k,j,i) = tend(k,j,i) - ukomp *                         &
110                                         ( v(k,j,i) - v(k,j,i-1) ) * ddx
111                ELSE
112                   tend(k,j,i) = tend(k,j,i) - ukomp *                         &
113                                         ( v(k,j,i+1) - v(k,j,i) ) * ddx
114                ENDIF
115!
116!--             y-direction
117                vkomp = v(k,j,i) - v_gtrans
118                IF ( vkomp > 0.0_wp )  THEN
119                   tend(k,j,i) = tend(k,j,i) - vkomp *                         &
120                                         ( v(k,j,i) - v(k,j-1,i) ) * ddy
121                ELSE
122                   tend(k,j,i) = tend(k,j,i) - vkomp *                         &
123                                         ( v(k,j+1,i) - v(k,j,i) ) * ddy
124                ENDIF
125!
126!--             z-direction
127                wkomp = 0.25_wp * ( w(k,j,i)  + w(k-1,j,i) +                   &
128                                 w(k,j-1,i) + w(k-1,j-1,i) )
129                IF ( wkomp > 0.0_wp )  THEN
130                   tend(k,j,i) = tend(k,j,i) - wkomp *                         &
131                                         ( v(k,j,i) - v(k-1,j,i) ) * ddzu(k)
132                ELSE
133                   tend(k,j,i) = tend(k,j,i) - wkomp *                         &
134                                         ( v(k+1,j,i) - v(k,j,i) ) * ddzu(k+1)
135                ENDIF
136
137             ENDDO
138          ENDDO
139       ENDDO
140
141    END SUBROUTINE advec_v_up
142
143
144!------------------------------------------------------------------------------!
145! Description:
146! ------------
147!> Call for grid point i,j
148!------------------------------------------------------------------------------!
149    SUBROUTINE advec_v_up_ij( i, j )
150
151       USE arrays_3d,                                                          &
152           ONLY:  ddzu, tend, u, v, w
153
154       USE control_parameters,                                                 &
155           ONLY:  u_gtrans, v_gtrans
156
157       USE grid_variables,                                                     &
158           ONLY:  ddx, ddy
159
160       USE indices,                                                            &
161           ONLY:  nzb_v_inner, nzt
162
163       USE kinds
164
165
166       IMPLICIT NONE
167
168       INTEGER(iwp) ::  i !<
169       INTEGER(iwp) ::  j !<
170       INTEGER(iwp) ::  k !<
171
172       REAL(wp) ::  ukomp !<
173       REAL(wp) ::  vkomp !<
174       REAL(wp) ::  wkomp !<
175
176
177       DO  k = nzb_v_inner(j,i)+1, nzt
178!
179!--       x-direction
180          ukomp = 0.25_wp * ( u(k,j,i) + u(k,j-1,i) + u(k,j,i+1) + u(k,j-1,i+1) &
181                         ) - u_gtrans
182          IF ( ukomp > 0.0_wp )  THEN
183             tend(k,j,i) = tend(k,j,i) - ukomp *                               &
184                                         ( v(k,j,i) - v(k,j,i-1) ) * ddx
185          ELSE
186             tend(k,j,i) = tend(k,j,i) - ukomp *                               &
187                                         ( v(k,j,i+1) - v(k,j,i) ) * ddx
188          ENDIF
189!
190!--       y-direction
191          vkomp = v(k,j,i) - v_gtrans
192          IF ( vkomp > 0.0_wp )  THEN
193             tend(k,j,i) = tend(k,j,i) - vkomp *                               &
194                                         ( v(k,j,i) - v(k,j-1,i) ) * ddy
195          ELSE
196             tend(k,j,i) = tend(k,j,i) - vkomp *                               &
197                                         ( v(k,j+1,i) - v(k,j,i) ) * ddy
198          ENDIF
199!
200!--       z-direction
201          wkomp = 0.25_wp * ( w(k,j,i) + w(k-1,j,i) + w(k,j-1,i) + w(k-1,j-1,i) )
202          IF ( wkomp > 0.0_wp )  THEN
203             tend(k,j,i) = tend(k,j,i) - wkomp *                               &
204                                         ( v(k,j,i) - v(k-1,j,i) ) * ddzu(k)
205          ELSE
206             tend(k,j,i) = tend(k,j,i) - wkomp *                               &
207                                         ( v(k+1,j,i) - v(k,j,i) ) * ddzu(k+1)
208          ENDIF
209
210       ENDDO
211
212    END SUBROUTINE advec_v_up_ij
213
214 END MODULE advec_v_up_mod
Note: See TracBrowser for help on using the repository browser.