source: palm/trunk/SOURCE/coriolis.f90 @ 1

Last change on this file since 1 was 1, checked in by raasch, 17 years ago

Initial repository layout and content

File size: 6.1 KB
Line 
1 MODULE coriolis_mod
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Log: coriolis.f90,v $
11! Revision 1.12  2006/02/23 10:08:57  raasch
12! nzb_2d replaced by nzb_u/v/w_inner
13!
14! Revision 1.11  2005/06/29 10:30:17  steinfeld
15! A potential dependency of the geostrophic wind components ug and vg on height
16! is considered
17!
18! Revision 1.10  2005/03/26 20:02:04  raasch
19! Extension of horizontal loop upper bounds for non-cyclic boundary conditions
20!
21! Revision 1.9  2004/01/30 10:17:34  raasch
22! Scalar lower k index nzb replaced by 2d-array nzb_2d
23!
24! Revision 1.8  2003/03/12 16:23:47  raasch
25! Full code replaced in the call for all gridpoints instead of calling the
26! _ij version (required by NEC, because otherwise no vectorization)
27!
28! Revision 1.7  2002/12/19 14:08:14  raasch
29! STOP statement replaced by call of subroutine local_stop
30!
31! Revision 1.6  2002/06/11 12:51:09  raasch
32! Former subroutine changed to a module which allows to be called for all grid
33! points of a single vertical column with index i,j or for all grid points by
34! using function overloading.
35!
36! Revision 1.5  2001/03/30 06:59:09  raasch
37! Translation of remaining German identifiers (variables, subroutines, etc.)
38!
39! Revision 1.4  2001/01/22 05:44:16  raasch
40! Module test_variables removed
41!
42! Revision 1.3  2000/12/20 09:50:27  letzel
43! All comments translated into English.
44!
45! Revision 1.2  1998/07/06 12:09:54  raasch
46! + USE test_variables
47!
48! Revision 1.1  1997/08/29 08:57:38  raasch
49! Initial revision
50!
51!
52! Description:
53! ------------
54! Computation of all Coriolis terms in the equations of motion.
55!------------------------------------------------------------------------------!
56
57    PRIVATE
58    PUBLIC coriolis
59
60    INTERFACE coriolis
61       MODULE PROCEDURE coriolis
62       MODULE PROCEDURE coriolis_ij
63    END INTERFACE coriolis
64
65 CONTAINS
66
67
68!------------------------------------------------------------------------------!
69! Call for all grid points
70!------------------------------------------------------------------------------!
71    SUBROUTINE coriolis( component )
72
73       USE arrays_3d
74       USE control_parameters
75       USE indices
76       USE pegrid
77
78       IMPLICIT NONE
79
80       INTEGER ::  component, i, j, k
81
82
83!
84!--    Compute Coriolis terms for the three velocity components
85       SELECT CASE ( component )
86
87!
88!--       u-component
89          CASE ( 1 )
90             DO  i = nxl, nxr+uxrp
91                DO  j = nys, nyn
92                   DO  k = nzb_u_inner(j,i)+1, nzt
93                      tend(k,j,i) = tend(k,j,i) + f  *    ( 0.25 *            &
94                                   ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +   &
95                                     v(k,j+1,i) ) - vg(k) )                   &
96                                             - fs *    ( 0.25 *               &
97                                   ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + &
98                                     w(k,j,i)   ) &
99                                                          )
100                   ENDDO
101                ENDDO
102             ENDDO
103
104!
105!--       v-component
106          CASE ( 2 )
107             DO  i = nxl, nxr
108                DO  j = nys, nyn+vynp
109                   DO  k = nzb_v_inner(j,i)+1, nzt
110                      tend(k,j,i) = tend(k,j,i) - f *     ( 0.25 *          &
111                                   ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + &
112                                     u(k,j,i+1) ) - ug(k) )
113                   ENDDO
114                ENDDO
115             ENDDO
116
117!
118!--       w-component
119          CASE ( 3 )
120             DO  i = nxl, nxr
121                DO  j = nys, nyn
122                   DO  k = nzb_w_inner(j,i)+1, nzt
123                      tend(k,j,i) = tend(k,j,i) + fs * 0.25 *             &
124                                   ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + &
125                                     u(k+1,j,i+1) )
126                   ENDDO
127                ENDDO
128             ENDDO
129
130          CASE DEFAULT
131
132             IF ( myid == 0 )  PRINT*,'+++ coriolis:  wrong component: ', &
133                                      component
134             CALL local_stop
135
136       END SELECT
137
138    END SUBROUTINE coriolis
139
140
141!------------------------------------------------------------------------------!
142! Call for grid point i,j
143!------------------------------------------------------------------------------!
144    SUBROUTINE coriolis_ij( i, j, component )
145
146       USE arrays_3d
147       USE control_parameters
148       USE indices
149       USE pegrid
150
151       IMPLICIT NONE
152
153       INTEGER ::  component, i, j, k
154
155!
156!--    Compute Coriolis terms for the three velocity components
157       SELECT CASE ( component )
158
159!
160!--       u-component
161          CASE ( 1 )
162             DO  k = nzb_u_inner(j,i)+1, nzt
163                tend(k,j,i) = tend(k,j,i) + f  *    ( 0.25 *               &
164                                ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +   &
165                                  v(k,j+1,i) ) - vg(k) )                   &
166                                          - fs *    ( 0.25 *               &
167                                ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + &
168                                  w(k,j,i)   ) &
169                                                    )
170             ENDDO
171
172!
173!--       v-component
174          CASE ( 2 )
175             DO  k = nzb_v_inner(j,i)+1, nzt
176                tend(k,j,i) = tend(k,j,i) - f *     ( 0.25 *             &
177                                ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + &
178                                  u(k,j,i+1) ) - ug(k) )
179             ENDDO
180
181!
182!--       w-component
183          CASE ( 3 )
184             DO  k = nzb_w_inner(j,i)+1, nzt
185                tend(k,j,i) = tend(k,j,i) + fs * 0.25 * &
186                                ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + &
187                                  u(k+1,j,i+1) )
188             ENDDO
189
190          CASE DEFAULT
191
192             IF ( myid == 0 )  PRINT*,'+++ coriolis:  wrong component: ', &
193                                      component
194             CALL local_stop
195
196       END SELECT
197
198    END SUBROUTINE coriolis_ij
199
200 END MODULE coriolis_mod
Note: See TracBrowser for help on using the repository browser.