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

Last change on this file since 3182 was 3182, checked in by suehring, 6 years ago

New Inifor features: grid stretching, improved command-interface, support start dates in different formats in both YYYYMMDD and YYYYMMDDHH, Ability to manually control input file prefixes (--radiation-prefix, --soil-preifx, --flow-prefix, --soilmoisture-prefix) for compatiblity with DWD forcast naming scheme; GNU-style short and long option; Prepared output of large-scale forcing profiles (no computation yet); Added preprocessor flag netcdf4 to switch output format between netCDF 3 and 4; Updated netCDF variable names and attributes to comply with PIDS v1.9; Inifor bugfixes: Improved compatibility with older Intel Intel compilers by avoiding implicit array allocation; Added origin_lon/_lat values and correct reference time in dynamic driver global attributes; corresponding PALM changes: adjustments to revised Inifor; variables names in dynamic driver adjusted; enable geostrophic forcing also in offline nested mode; variable names in LES-LES and COSMO offline nesting changed; lateral boundary flags for nesting, in- and outflow conditions renamed

  • Property svn:keywords set to Id
File size: 9.4 KB
Line 
1!> @file coriolis.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
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.
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!
17! Copyright 1997-2018 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22! Remove masking of geostrophic wind forcing in offline nesting case
23!
24! Former revisions:
25! -----------------
26! $Id: coriolis.f90 3182 2018-07-27 13:36:03Z suehring $
27! Corrected "Former revisions" section
28!
29! 2696 2017-12-14 17:12:51Z kanani
30! Change in file header (GPL part)
31! Forcing implemented, preliminary (MS)
32!
33! 2233 2017-05-30 18:08:54Z suehring
34!
35! 2232 2017-05-30 17:47:52Z suehring
36! Adjustments to new topography concept
37!
38! 2118 2017-01-17 16:38:49Z raasch
39! OpenACC version of subroutine removed
40!
41! 2000 2016-08-20 18:09:15Z knoop
42! Forced header and separation lines into 80 columns
43!
44! 1873 2016-04-18 14:50:06Z maronga
45! Module renamed (removed _mod)
46!
47!
48! 1850 2016-04-08 13:29:27Z maronga
49! Module renamed
50!
51! 1682 2015-10-07 23:56:08Z knoop
52! Code annotations made doxygen readable
53!
54! 1353 2014-04-08 15:21:23Z heinze
55! REAL constants provided with KIND-attribute
56!
57! 1320 2014-03-20 08:40:49Z raasch
58! ONLY-attribute added to USE-statements,
59! kind-parameters added to all INTEGER and REAL declaration statements,
60! kinds are defined in new module kinds,
61! revision history before 2012 removed,
62! comment fields (!:) to be used for variable explanations added to
63! all variable declaration statements
64!
65! 1257 2013-11-08 15:18:40Z raasch
66! openacc loop and loop vector clauses removed
67!
68! 1128 2013-04-12 06:19:32Z raasch
69! loop index bounds in accelerator version replaced by i_left, i_right, j_south,
70! j_north
71!
72! 1036 2012-10-22 13:43:42Z raasch
73! code put under GPL (PALM 3.9)
74!
75! 1015 2012-09-27 09:23:24Z raasch
76! accelerator version (*_acc) added
77!
78! Revision 1.1  1997/08/29 08:57:38  raasch
79! Initial revision
80!
81!
82! Description:
83! ------------
84!> Computation of all Coriolis terms in the equations of motion.
85!------------------------------------------------------------------------------!
86 MODULE coriolis_mod
87 
88
89    PRIVATE
90    PUBLIC coriolis
91
92    INTERFACE coriolis
93       MODULE PROCEDURE coriolis
94       MODULE PROCEDURE coriolis_ij
95    END INTERFACE coriolis
96
97 CONTAINS
98
99
100!------------------------------------------------------------------------------!
101! Description:
102! ------------
103!> Call for all grid points
104!------------------------------------------------------------------------------!
105    SUBROUTINE coriolis( component )
106
107       USE arrays_3d,                                                          &
108           ONLY:  tend, u, ug, v, vg, w 
109           
110       USE control_parameters,                                                 &
111           ONLY:  f, fs, message_string, nesting_offline
112           
113       USE indices,                                                            &
114           ONLY:  nxl, nxlu, nxr, nyn, nys, nysv, nzb, nzt, wall_flags_0
115                   
116       USE kinds
117
118       IMPLICIT NONE
119
120       INTEGER(iwp) ::  component  !<
121       INTEGER(iwp) ::  i          !< running index x direction
122       INTEGER(iwp) ::  j          !< running index y direction
123       INTEGER(iwp) ::  k          !< running index z direction
124
125       REAL(wp)     ::  flag           !< flag to mask topography
126
127!
128!--    Compute Coriolis terms for the three velocity components
129       SELECT CASE ( component )
130
131!
132!--       u-component
133          CASE ( 1 )
134             DO  i = nxlu, nxr
135                DO  j = nys, nyn
136                   DO  k = nzb+1, nzt
137!
138!--                   Predetermine flag to mask topography
139                      flag = MERGE( 1.0_wp, 0.0_wp,                            &
140                                    BTEST( wall_flags_0(k,j,i), 1 ) )
141
142                      tend(k,j,i) = tend(k,j,i) + f  *    ( 0.25_wp *          &
143                                   ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +    &
144                                     v(k,j+1,i) ) - vg(k) ) * flag             &
145                                                - fs *    ( 0.25_wp *          &
146                                   ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) +  &
147                                     w(k,j,i)   )                              &
148                                                          ) * flag
149                   ENDDO
150                ENDDO
151             ENDDO
152
153!
154!--       v-component
155          CASE ( 2 )
156             DO  i = nxl, nxr
157                DO  j = nysv, nyn
158                   DO  k = nzb+1, nzt
159!
160!--                   Predetermine flag to mask topography
161                      flag = MERGE( 1.0_wp, 0.0_wp,                            &
162                                    BTEST( wall_flags_0(k,j,i), 2 ) )
163
164                      tend(k,j,i) = tend(k,j,i) - f *     ( 0.25_wp *          &
165                                   ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) +    &
166                                     u(k,j,i+1) ) - ug(k) ) * flag
167                   ENDDO
168                ENDDO
169             ENDDO
170
171!
172!--       w-component
173          CASE ( 3 )
174             DO  i = nxl, nxr
175                DO  j = nys, nyn
176                   DO  k = nzb+1, nzt
177!
178!--                   Predetermine flag to mask topography
179                      flag = MERGE( 1.0_wp, 0.0_wp,                            &
180                                    BTEST( wall_flags_0(k,j,i), 3 ) )
181
182                      tend(k,j,i) = tend(k,j,i) + fs * 0.25_wp *               &
183                                   ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) +      &
184                                     u(k+1,j,i+1) ) * flag
185                   ENDDO
186                ENDDO
187             ENDDO
188
189          CASE DEFAULT
190
191             WRITE( message_string, * ) ' wrong component: ', component
192             CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 )
193
194       END SELECT
195
196    END SUBROUTINE coriolis
197
198
199!------------------------------------------------------------------------------!
200! Description:
201! ------------
202!> Call for grid point i,j
203!------------------------------------------------------------------------------!
204    SUBROUTINE coriolis_ij( i, j, component )
205
206       USE arrays_3d,                                                          &
207           ONLY:  tend, u, ug, v, vg, w 
208           
209       USE control_parameters,                                                 &
210           ONLY:  f, fs, message_string, nesting_offline
211           
212       USE indices,                                                            &
213           ONLY:  nzb, nzt, wall_flags_0
214           
215       USE kinds
216       
217       IMPLICIT NONE
218
219       INTEGER(iwp) ::  component  !<
220       INTEGER(iwp) ::  i          !< running index x direction
221       INTEGER(iwp) ::  j          !< running index y direction
222       INTEGER(iwp) ::  k          !< running index z direction
223
224       REAL(wp)     ::  flag       !< flag to mask topography
225
226!
227!--    Compute Coriolis terms for the three velocity components
228       SELECT CASE ( component )
229
230!
231!--       u-component
232          CASE ( 1 )
233             DO  k = nzb+1, nzt
234!
235!--             Predetermine flag to mask topography
236                flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 1 ) )
237
238                tend(k,j,i) = tend(k,j,i) + f  *     ( 0.25_wp *               &
239                                ( v(k,j,i-1) + v(k,j,i) + v(k,j+1,i-1) +       &
240                                  v(k,j+1,i) ) - vg(k)                         &
241                                                     ) * flag                  &
242                                          - fs *     ( 0.25_wp *               &
243                                ( w(k-1,j,i-1) + w(k-1,j,i) + w(k,j,i-1) +     &
244                                  w(k,j,i)   )       ) * flag
245             ENDDO
246
247!
248!--       v-component
249          CASE ( 2 )
250             DO  k = nzb+1, nzt
251!
252!--             Predetermine flag to mask topography
253                flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 2 ) )
254
255                tend(k,j,i) = tend(k,j,i) - f *        ( 0.25_wp *             &
256                                ( u(k,j-1,i) + u(k,j,i) + u(k,j-1,i+1) +       &
257                                  u(k,j,i+1) ) - ug(k) ) * flag
258             ENDDO
259
260!
261!--       w-component
262          CASE ( 3 )
263             DO  k = nzb+1, nzt
264!
265!--             Predetermine flag to mask topography
266                flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 3 ) )
267
268                tend(k,j,i) = tend(k,j,i) + fs * 0.25_wp *                     &
269                                ( u(k,j,i) + u(k+1,j,i) + u(k,j,i+1) +         &
270                                  u(k+1,j,i+1) ) * flag
271             ENDDO
272
273          CASE DEFAULT
274
275             WRITE( message_string, * ) ' wrong component: ', component
276             CALL message( 'coriolis', 'PA0173', 1, 2, 0, 6, 0 )
277
278       END SELECT
279
280    END SUBROUTINE coriolis_ij
281
282 END MODULE coriolis_mod
Note: See TracBrowser for help on using the repository browser.