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

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

preliminary update of bugfixes and extensions for non-cyclic BCs

  • Property svn:keywords set to Id
File size: 5.1 KB
Line 
1 MODULE coriolis_mod
2
3!------------------------------------------------------------------------------!
4! Actual revisions:
5! -----------------
6! loops for u and v are starting from index nxlu, nysv, respectively (needed
7! for non-cyclic boundary conditions)
8!
9! Former revisions:
10! -----------------
11! $Id: coriolis.f90 106 2007-08-16 14:30:26Z raasch $
12!
13! 75 2007-03-22 09:54:05Z raasch
14! uxrp, vynp eliminated
15!
16! RCS Log replace by Id keyword, revision history cleaned up
17!
18! Revision 1.12  2006/02/23 10:08:57  raasch
19! nzb_2d replaced by nzb_u/v/w_inner
20!
21! Revision 1.1  1997/08/29 08:57:38  raasch
22! Initial revision
23!
24!
25! Description:
26! ------------
27! Computation of all Coriolis terms in the equations of motion.
28!------------------------------------------------------------------------------!
29
30    PRIVATE
31    PUBLIC coriolis
32
33    INTERFACE coriolis
34       MODULE PROCEDURE coriolis
35       MODULE PROCEDURE coriolis_ij
36    END INTERFACE coriolis
37
38 CONTAINS
39
40
41!------------------------------------------------------------------------------!
42! Call for all grid points
43!------------------------------------------------------------------------------!
44    SUBROUTINE coriolis( component )
45
46       USE arrays_3d
47       USE control_parameters
48       USE indices
49       USE pegrid
50
51       IMPLICIT NONE
52
53       INTEGER ::  component, i, j, k
54
55
56!
57!--    Compute Coriolis terms for the three velocity components
58       SELECT CASE ( component )
59
60!
61!--       u-component
62          CASE ( 1 )
63             DO  i = nxlu, nxr
64                DO  j = nys, nyn
65                   DO  k = nzb_u_inner(j,i)+1, nzt
66                      tend(k,j,i) = tend(k,j,i) + f  *    ( 0.25 *            &
67                                   ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +   &
68                                     v(k,j+1,i) ) - vg(k) )                   &
69                                             - fs *    ( 0.25 *               &
70                                   ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + &
71                                     w(k,j,i)   ) &
72                                                          )
73                   ENDDO
74                ENDDO
75             ENDDO
76
77!
78!--       v-component
79          CASE ( 2 )
80             DO  i = nxl, nxr
81                DO  j = nysv, nyn
82                   DO  k = nzb_v_inner(j,i)+1, nzt
83                      tend(k,j,i) = tend(k,j,i) - f *     ( 0.25 *          &
84                                   ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + &
85                                     u(k,j,i+1) ) - ug(k) )
86                   ENDDO
87                ENDDO
88             ENDDO
89
90!
91!--       w-component
92          CASE ( 3 )
93             DO  i = nxl, nxr
94                DO  j = nys, nyn
95                   DO  k = nzb_w_inner(j,i)+1, nzt
96                      tend(k,j,i) = tend(k,j,i) + fs * 0.25 *             &
97                                   ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + &
98                                     u(k+1,j,i+1) )
99                   ENDDO
100                ENDDO
101             ENDDO
102
103          CASE DEFAULT
104
105             IF ( myid == 0 )  PRINT*,'+++ coriolis:  wrong component: ', &
106                                      component
107             CALL local_stop
108
109       END SELECT
110
111    END SUBROUTINE coriolis
112
113
114!------------------------------------------------------------------------------!
115! Call for grid point i,j
116!------------------------------------------------------------------------------!
117    SUBROUTINE coriolis_ij( i, j, component )
118
119       USE arrays_3d
120       USE control_parameters
121       USE indices
122       USE pegrid
123
124       IMPLICIT NONE
125
126       INTEGER ::  component, i, j, k
127
128!
129!--    Compute Coriolis terms for the three velocity components
130       SELECT CASE ( component )
131
132!
133!--       u-component
134          CASE ( 1 )
135             DO  k = nzb_u_inner(j,i)+1, nzt
136                tend(k,j,i) = tend(k,j,i) + f  *    ( 0.25 *               &
137                                ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +   &
138                                  v(k,j+1,i) ) - vg(k) )                   &
139                                          - fs *    ( 0.25 *               &
140                                ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) + &
141                                  w(k,j,i)   ) &
142                                                    )
143             ENDDO
144
145!
146!--       v-component
147          CASE ( 2 )
148             DO  k = nzb_v_inner(j,i)+1, nzt
149                tend(k,j,i) = tend(k,j,i) - f *     ( 0.25 *             &
150                                ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) + &
151                                  u(k,j,i+1) ) - ug(k) )
152             ENDDO
153
154!
155!--       w-component
156          CASE ( 3 )
157             DO  k = nzb_w_inner(j,i)+1, nzt
158                tend(k,j,i) = tend(k,j,i) + fs * 0.25 * &
159                                ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) + &
160                                  u(k+1,j,i+1) )
161             ENDDO
162
163          CASE DEFAULT
164
165             IF ( myid == 0 )  PRINT*,'+++ coriolis:  wrong component: ', &
166                                      component
167             CALL local_stop
168
169       END SELECT
170
171    END SUBROUTINE coriolis_ij
172
173 END MODULE coriolis_mod
Note: See TracBrowser for help on using the repository browser.