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

Last change on this file since 2803 was 2718, checked in by maronga, 7 years ago

deleting of deprecated files; headers updated where needed

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