[1682] | 1 | !> @file boundary_conds.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[1036] | 3 | ! This file is part of PALM. |
---|
| 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 | ! |
---|
[2101] | 17 | ! Copyright 1997-2017 Leibniz Universitaet Hannover |
---|
[2000] | 18 | !------------------------------------------------------------------------------! |
---|
[1036] | 19 | ! |
---|
[484] | 20 | ! Current revisions: |
---|
[1] | 21 | ! ----------------- |
---|
[1933] | 22 | ! |
---|
[2119] | 23 | ! |
---|
[1321] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: boundary_conds.f90 2119 2017-01-17 16:51:50Z lvollmer $ |
---|
| 27 | ! |
---|
[2119] | 28 | ! 2118 2017-01-17 16:38:49Z raasch |
---|
| 29 | ! OpenACC directives removed |
---|
| 30 | ! |
---|
[2001] | 31 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
| 32 | ! Forced header and separation lines into 80 columns |
---|
| 33 | ! |
---|
[1993] | 34 | ! 1992 2016-08-12 15:14:59Z suehring |
---|
| 35 | ! Adjustments for top boundary condition for passive scalar |
---|
| 36 | ! |
---|
[1961] | 37 | ! 1960 2016-07-12 16:34:24Z suehring |
---|
| 38 | ! Treat humidity and passive scalar separately |
---|
| 39 | ! |
---|
[1933] | 40 | ! 1823 2016-04-07 08:57:52Z hoffmann |
---|
| 41 | ! Initial version of purely vertical nesting introduced. |
---|
| 42 | ! |
---|
[1823] | 43 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
| 44 | ! icloud_scheme removed. microphyisics_seifert added. |
---|
| 45 | ! |
---|
[1765] | 46 | ! 1764 2016-02-28 12:45:19Z raasch |
---|
| 47 | ! index bug for u_p at left outflow removed |
---|
| 48 | ! |
---|
[1763] | 49 | ! 1762 2016-02-25 12:31:13Z hellstea |
---|
| 50 | ! Introduction of nested domain feature |
---|
| 51 | ! |
---|
[1744] | 52 | ! 1742 2016-01-13 09:50:06Z raasch |
---|
| 53 | ! bugfix for outflow Neumann boundary conditions at bottom and top |
---|
| 54 | ! |
---|
[1718] | 55 | ! 1717 2015-11-11 15:09:47Z raasch |
---|
| 56 | ! Bugfix: index error in outflow conditions for left boundary |
---|
| 57 | ! |
---|
[1683] | 58 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 59 | ! Code annotations made doxygen readable |
---|
| 60 | ! |
---|
[1717] | 61 | ! 1410 2014-05-23 12:16:18Z suehring |
---|
[1463] | 62 | ! Bugfix: set dirichlet boundary condition for passive_scalar at model domain |
---|
| 63 | ! top |
---|
| 64 | ! |
---|
[1410] | 65 | ! 1399 2014-05-07 11:16:25Z heinze |
---|
| 66 | ! Bugfix: set inflow boundary conditions also if no humidity or passive_scalar |
---|
| 67 | ! is used. |
---|
| 68 | ! |
---|
[1399] | 69 | ! 1398 2014-05-07 11:15:00Z heinze |
---|
| 70 | ! Dirichlet-condition at the top for u and v changed to u_init and v_init also |
---|
| 71 | ! for large_scale_forcing |
---|
| 72 | ! |
---|
[1381] | 73 | ! 1380 2014-04-28 12:40:45Z heinze |
---|
| 74 | ! Adjust Dirichlet-condition at the top for pt in case of nudging |
---|
| 75 | ! |
---|
[1362] | 76 | ! 1361 2014-04-16 15:17:48Z hoffmann |
---|
| 77 | ! Bottom and top boundary conditions of rain water content (qr) and |
---|
| 78 | ! rain drop concentration (nr) changed to Dirichlet |
---|
| 79 | ! |
---|
[1354] | 80 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
| 81 | ! REAL constants provided with KIND-attribute |
---|
| 82 | ! |
---|
[1321] | 83 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 84 | ! ONLY-attribute added to USE-statements, |
---|
| 85 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 86 | ! kinds are defined in new module kinds, |
---|
| 87 | ! revision history before 2012 removed, |
---|
| 88 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 89 | ! all variable declaration statements |
---|
[1160] | 90 | ! |
---|
[1258] | 91 | ! 1257 2013-11-08 15:18:40Z raasch |
---|
| 92 | ! loop independent clauses added |
---|
| 93 | ! |
---|
[1242] | 94 | ! 1241 2013-10-30 11:36:58Z heinze |
---|
| 95 | ! Adjust ug and vg at each timestep in case of large_scale_forcing |
---|
| 96 | ! |
---|
[1160] | 97 | ! 1159 2013-05-21 11:58:22Z fricke |
---|
[1159] | 98 | ! Bugfix: Neumann boundary conditions for the velocity components at the |
---|
| 99 | ! outflow are in fact radiation boundary conditions using the maximum phase |
---|
| 100 | ! velocity that ensures numerical stability (CFL-condition). |
---|
| 101 | ! Hence, logical operator use_cmax is now used instead of bc_lr_dirneu/_neudir. |
---|
| 102 | ! Bugfix: In case of use_cmax at the outflow, u, v, w are replaced by |
---|
| 103 | ! u_p, v_p, w_p |
---|
[1116] | 104 | ! |
---|
| 105 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
| 106 | ! boundary conditions of two-moment cloud scheme are restricted to Neumann- |
---|
| 107 | ! boundary-conditions |
---|
| 108 | ! |
---|
[1114] | 109 | ! 1113 2013-03-10 02:48:14Z raasch |
---|
| 110 | ! GPU-porting |
---|
| 111 | ! dummy argument "range" removed |
---|
| 112 | ! Bugfix: wrong index in loops of radiation boundary condition |
---|
[1113] | 113 | ! |
---|
[1054] | 114 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
| 115 | ! boundary conditions for the two new prognostic equations (nr, qr) of the |
---|
| 116 | ! two-moment cloud scheme |
---|
| 117 | ! |
---|
[1037] | 118 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 119 | ! code put under GPL (PALM 3.9) |
---|
| 120 | ! |
---|
[997] | 121 | ! 996 2012-09-07 10:41:47Z raasch |
---|
| 122 | ! little reformatting |
---|
| 123 | ! |
---|
[979] | 124 | ! 978 2012-08-09 08:28:32Z fricke |
---|
| 125 | ! Neumann boudnary conditions are added at the inflow boundary for the SGS-TKE. |
---|
| 126 | ! Outflow boundary conditions for the velocity components can be set to Neumann |
---|
| 127 | ! conditions or to radiation conditions with a horizontal averaged phase |
---|
| 128 | ! velocity. |
---|
| 129 | ! |
---|
[876] | 130 | ! 875 2012-04-02 15:35:15Z gryschka |
---|
| 131 | ! Bugfix in case of dirichlet inflow bc at the right or north boundary |
---|
| 132 | ! |
---|
[1] | 133 | ! Revision 1.1 1997/09/12 06:21:34 raasch |
---|
| 134 | ! Initial revision |
---|
| 135 | ! |
---|
| 136 | ! |
---|
| 137 | ! Description: |
---|
| 138 | ! ------------ |
---|
[1682] | 139 | !> Boundary conditions for the prognostic quantities. |
---|
| 140 | !> One additional bottom boundary condition is applied for the TKE (=(u*)**2) |
---|
| 141 | !> in prandtl_fluxes. The cyclic lateral boundary conditions are implicitly |
---|
| 142 | !> handled in routine exchange_horiz. Pressure boundary conditions are |
---|
| 143 | !> explicitly set in routines pres, poisfft, poismg and sor. |
---|
[1] | 144 | !------------------------------------------------------------------------------! |
---|
[1682] | 145 | SUBROUTINE boundary_conds |
---|
| 146 | |
---|
[1] | 147 | |
---|
[1320] | 148 | USE arrays_3d, & |
---|
| 149 | ONLY: c_u, c_u_m, c_u_m_l, c_v, c_v_m, c_v_m_l, c_w, c_w_m, c_w_m_l, & |
---|
[1960] | 150 | dzu, e_p, nr_p, pt, pt_p, q, q_p, qr_p, s, s_p, sa, sa_p, & |
---|
[1320] | 151 | u, ug, u_init, u_m_l, u_m_n, u_m_r, u_m_s, u_p, & |
---|
| 152 | v, vg, v_init, v_m_l, v_m_n, v_m_r, v_m_s, v_p, & |
---|
[1960] | 153 | w, w_p, w_m_l, w_m_n, w_m_r, w_m_s, pt_init |
---|
[1320] | 154 | |
---|
| 155 | USE control_parameters, & |
---|
[1960] | 156 | ONLY: bc_pt_t_val, bc_q_t_val, bc_s_t_val, constant_diffusion, & |
---|
[1320] | 157 | cloud_physics, dt_3d, humidity, & |
---|
[1960] | 158 | ibc_pt_b, ibc_pt_t, ibc_q_b, ibc_q_t, ibc_s_b, ibc_s_t, & |
---|
| 159 | ibc_sa_t, ibc_uv_b, ibc_uv_t, inflow_l, inflow_n, inflow_r, & |
---|
| 160 | inflow_s, intermediate_timestep_count, large_scale_forcing, & |
---|
[1822] | 161 | microphysics_seifert, nest_domain, nest_bound_l, nest_bound_s, & |
---|
| 162 | nudging, ocean, outflow_l, outflow_n, outflow_r, outflow_s, & |
---|
| 163 | passive_scalar, tsc, use_cmax |
---|
[1320] | 164 | |
---|
| 165 | USE grid_variables, & |
---|
| 166 | ONLY: ddx, ddy, dx, dy |
---|
| 167 | |
---|
| 168 | USE indices, & |
---|
| 169 | ONLY: nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, & |
---|
| 170 | nzb, nzb_s_inner, nzb_w_inner, nzt |
---|
| 171 | |
---|
| 172 | USE kinds |
---|
| 173 | |
---|
[1] | 174 | USE pegrid |
---|
| 175 | |
---|
[1933] | 176 | USE pmc_interface, & |
---|
| 177 | ONLY : nesting_mode |
---|
[1320] | 178 | |
---|
[1933] | 179 | |
---|
[1] | 180 | IMPLICIT NONE |
---|
| 181 | |
---|
[1682] | 182 | INTEGER(iwp) :: i !< |
---|
| 183 | INTEGER(iwp) :: j !< |
---|
| 184 | INTEGER(iwp) :: k !< |
---|
[1] | 185 | |
---|
[1682] | 186 | REAL(wp) :: c_max !< |
---|
| 187 | REAL(wp) :: denom !< |
---|
[1] | 188 | |
---|
[73] | 189 | |
---|
[1] | 190 | ! |
---|
[1113] | 191 | !-- Bottom boundary |
---|
| 192 | IF ( ibc_uv_b == 1 ) THEN |
---|
| 193 | u_p(nzb,:,:) = u_p(nzb+1,:,:) |
---|
| 194 | v_p(nzb,:,:) = v_p(nzb+1,:,:) |
---|
| 195 | ENDIF |
---|
| 196 | |
---|
| 197 | DO i = nxlg, nxrg |
---|
| 198 | DO j = nysg, nyng |
---|
[1353] | 199 | w_p(nzb_w_inner(j,i),j,i) = 0.0_wp |
---|
[1113] | 200 | ENDDO |
---|
| 201 | ENDDO |
---|
| 202 | |
---|
| 203 | ! |
---|
[1762] | 204 | !-- Top boundary. A nested domain ( ibc_uv_t = 3 ) does not require settings. |
---|
[1113] | 205 | IF ( ibc_uv_t == 0 ) THEN |
---|
| 206 | u_p(nzt+1,:,:) = u_init(nzt+1) |
---|
| 207 | v_p(nzt+1,:,:) = v_init(nzt+1) |
---|
[1762] | 208 | ELSEIF ( ibc_uv_t == 1 ) THEN |
---|
[1113] | 209 | u_p(nzt+1,:,:) = u_p(nzt,:,:) |
---|
| 210 | v_p(nzt+1,:,:) = v_p(nzt,:,:) |
---|
| 211 | ENDIF |
---|
| 212 | |
---|
[1762] | 213 | IF ( .NOT. nest_domain ) THEN |
---|
| 214 | w_p(nzt:nzt+1,:,:) = 0.0_wp ! nzt is not a prognostic level (but cf. pres) |
---|
| 215 | ENDIF |
---|
| 216 | |
---|
[1113] | 217 | ! |
---|
| 218 | !-- Temperature at bottom boundary. |
---|
| 219 | !-- In case of coupled runs (ibc_pt_b = 2) the temperature is given by |
---|
| 220 | !-- the sea surface temperature of the coupled ocean model. |
---|
| 221 | IF ( ibc_pt_b == 0 ) THEN |
---|
[667] | 222 | DO i = nxlg, nxrg |
---|
| 223 | DO j = nysg, nyng |
---|
[1113] | 224 | pt_p(nzb_s_inner(j,i),j,i) = pt(nzb_s_inner(j,i),j,i) |
---|
[1] | 225 | ENDDO |
---|
| 226 | ENDDO |
---|
[1113] | 227 | ELSEIF ( ibc_pt_b == 1 ) THEN |
---|
| 228 | DO i = nxlg, nxrg |
---|
| 229 | DO j = nysg, nyng |
---|
| 230 | pt_p(nzb_s_inner(j,i),j,i) = pt_p(nzb_s_inner(j,i)+1,j,i) |
---|
| 231 | ENDDO |
---|
| 232 | ENDDO |
---|
| 233 | ENDIF |
---|
[1] | 234 | |
---|
| 235 | ! |
---|
[1113] | 236 | !-- Temperature at top boundary |
---|
| 237 | IF ( ibc_pt_t == 0 ) THEN |
---|
| 238 | pt_p(nzt+1,:,:) = pt(nzt+1,:,:) |
---|
[1380] | 239 | ! |
---|
| 240 | !-- In case of nudging adjust top boundary to pt which is |
---|
| 241 | !-- read in from NUDGING-DATA |
---|
| 242 | IF ( nudging ) THEN |
---|
| 243 | pt_p(nzt+1,:,:) = pt_init(nzt+1) |
---|
| 244 | ENDIF |
---|
[1113] | 245 | ELSEIF ( ibc_pt_t == 1 ) THEN |
---|
| 246 | pt_p(nzt+1,:,:) = pt_p(nzt,:,:) |
---|
| 247 | ELSEIF ( ibc_pt_t == 2 ) THEN |
---|
[1992] | 248 | pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1) |
---|
[1113] | 249 | ENDIF |
---|
[1] | 250 | |
---|
| 251 | ! |
---|
[1113] | 252 | !-- Boundary conditions for TKE |
---|
| 253 | !-- Generally Neumann conditions with de/dz=0 are assumed |
---|
| 254 | IF ( .NOT. constant_diffusion ) THEN |
---|
| 255 | DO i = nxlg, nxrg |
---|
| 256 | DO j = nysg, nyng |
---|
| 257 | e_p(nzb_s_inner(j,i),j,i) = e_p(nzb_s_inner(j,i)+1,j,i) |
---|
[73] | 258 | ENDDO |
---|
[1113] | 259 | ENDDO |
---|
[1762] | 260 | IF ( .NOT. nest_domain ) THEN |
---|
| 261 | e_p(nzt+1,:,:) = e_p(nzt,:,:) |
---|
| 262 | ENDIF |
---|
[1113] | 263 | ENDIF |
---|
| 264 | |
---|
| 265 | ! |
---|
| 266 | !-- Boundary conditions for salinity |
---|
| 267 | IF ( ocean ) THEN |
---|
| 268 | ! |
---|
| 269 | !-- Bottom boundary: Neumann condition because salinity flux is always |
---|
| 270 | !-- given |
---|
| 271 | DO i = nxlg, nxrg |
---|
| 272 | DO j = nysg, nyng |
---|
| 273 | sa_p(nzb_s_inner(j,i),j,i) = sa_p(nzb_s_inner(j,i)+1,j,i) |
---|
[1] | 274 | ENDDO |
---|
[1113] | 275 | ENDDO |
---|
[1] | 276 | |
---|
| 277 | ! |
---|
[1113] | 278 | !-- Top boundary: Dirichlet or Neumann |
---|
| 279 | IF ( ibc_sa_t == 0 ) THEN |
---|
| 280 | sa_p(nzt+1,:,:) = sa(nzt+1,:,:) |
---|
| 281 | ELSEIF ( ibc_sa_t == 1 ) THEN |
---|
| 282 | sa_p(nzt+1,:,:) = sa_p(nzt,:,:) |
---|
[1] | 283 | ENDIF |
---|
| 284 | |
---|
[1113] | 285 | ENDIF |
---|
| 286 | |
---|
[1] | 287 | ! |
---|
[1960] | 288 | !-- Boundary conditions for total water content, |
---|
[1113] | 289 | !-- bottom and top boundary (see also temperature) |
---|
[1960] | 290 | IF ( humidity ) THEN |
---|
[1113] | 291 | ! |
---|
| 292 | !-- Surface conditions for constant_humidity_flux |
---|
| 293 | IF ( ibc_q_b == 0 ) THEN |
---|
[667] | 294 | DO i = nxlg, nxrg |
---|
| 295 | DO j = nysg, nyng |
---|
[1113] | 296 | q_p(nzb_s_inner(j,i),j,i) = q(nzb_s_inner(j,i),j,i) |
---|
[1] | 297 | ENDDO |
---|
| 298 | ENDDO |
---|
[1113] | 299 | ELSE |
---|
[667] | 300 | DO i = nxlg, nxrg |
---|
| 301 | DO j = nysg, nyng |
---|
[1113] | 302 | q_p(nzb_s_inner(j,i),j,i) = q_p(nzb_s_inner(j,i)+1,j,i) |
---|
[95] | 303 | ENDDO |
---|
| 304 | ENDDO |
---|
[1113] | 305 | ENDIF |
---|
[95] | 306 | ! |
---|
[1113] | 307 | !-- Top boundary |
---|
[1462] | 308 | IF ( ibc_q_t == 0 ) THEN |
---|
| 309 | q_p(nzt+1,:,:) = q(nzt+1,:,:) |
---|
| 310 | ELSEIF ( ibc_q_t == 1 ) THEN |
---|
[1992] | 311 | q_p(nzt+1,:,:) = q_p(nzt,:,:) + bc_q_t_val * dzu(nzt+1) |
---|
[1462] | 312 | ENDIF |
---|
[95] | 313 | |
---|
[1822] | 314 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
[1113] | 315 | ! |
---|
[1361] | 316 | !-- Surface conditions rain water (Dirichlet) |
---|
[1115] | 317 | DO i = nxlg, nxrg |
---|
| 318 | DO j = nysg, nyng |
---|
[1361] | 319 | qr_p(nzb_s_inner(j,i),j,i) = 0.0_wp |
---|
| 320 | nr_p(nzb_s_inner(j,i),j,i) = 0.0_wp |
---|
[73] | 321 | ENDDO |
---|
[1115] | 322 | ENDDO |
---|
[1] | 323 | ! |
---|
[1361] | 324 | !-- Top boundary condition for rain water (Dirichlet) |
---|
| 325 | qr_p(nzt+1,:,:) = 0.0_wp |
---|
| 326 | nr_p(nzt+1,:,:) = 0.0_wp |
---|
[1115] | 327 | |
---|
[1] | 328 | ENDIF |
---|
[1409] | 329 | ENDIF |
---|
[1] | 330 | ! |
---|
[1960] | 331 | !-- Boundary conditions for scalar, |
---|
| 332 | !-- bottom and top boundary (see also temperature) |
---|
| 333 | IF ( passive_scalar ) THEN |
---|
| 334 | ! |
---|
| 335 | !-- Surface conditions for constant_humidity_flux |
---|
| 336 | IF ( ibc_s_b == 0 ) THEN |
---|
| 337 | DO i = nxlg, nxrg |
---|
| 338 | DO j = nysg, nyng |
---|
| 339 | s_p(nzb_s_inner(j,i),j,i) = s(nzb_s_inner(j,i),j,i) |
---|
| 340 | ENDDO |
---|
| 341 | ENDDO |
---|
| 342 | ELSE |
---|
| 343 | DO i = nxlg, nxrg |
---|
| 344 | DO j = nysg, nyng |
---|
| 345 | s_p(nzb_s_inner(j,i),j,i) = s_p(nzb_s_inner(j,i)+1,j,i) |
---|
| 346 | ENDDO |
---|
| 347 | ENDDO |
---|
| 348 | ENDIF |
---|
| 349 | ! |
---|
[1992] | 350 | !-- Top boundary condition |
---|
| 351 | IF ( ibc_s_t == 0 ) THEN |
---|
[1960] | 352 | s_p(nzt+1,:,:) = s(nzt+1,:,:) |
---|
[1992] | 353 | ELSEIF ( ibc_s_t == 1 ) THEN |
---|
| 354 | s_p(nzt+1,:,:) = s_p(nzt,:,:) |
---|
| 355 | ELSEIF ( ibc_s_t == 2 ) THEN |
---|
| 356 | s_p(nzt+1,:,:) = s_p(nzt,:,:) + bc_s_t_val * dzu(nzt+1) |
---|
[1960] | 357 | ENDIF |
---|
| 358 | |
---|
| 359 | ENDIF |
---|
| 360 | ! |
---|
[1762] | 361 | !-- In case of inflow or nest boundary at the south boundary the boundary for v |
---|
| 362 | !-- is at nys and in case of inflow or nest boundary at the left boundary the |
---|
| 363 | !-- boundary for u is at nxl. Since in prognostic_equations (cache optimized |
---|
| 364 | !-- version) these levels are handled as a prognostic level, boundary values |
---|
| 365 | !-- have to be restored here. |
---|
[1409] | 366 | !-- For the SGS-TKE, Neumann boundary conditions are used at the inflow. |
---|
| 367 | IF ( inflow_s ) THEN |
---|
| 368 | v_p(:,nys,:) = v_p(:,nys-1,:) |
---|
| 369 | IF ( .NOT. constant_diffusion ) e_p(:,nys-1,:) = e_p(:,nys,:) |
---|
| 370 | ELSEIF ( inflow_n ) THEN |
---|
| 371 | IF ( .NOT. constant_diffusion ) e_p(:,nyn+1,:) = e_p(:,nyn,:) |
---|
| 372 | ELSEIF ( inflow_l ) THEN |
---|
| 373 | u_p(:,:,nxl) = u_p(:,:,nxl-1) |
---|
| 374 | IF ( .NOT. constant_diffusion ) e_p(:,:,nxl-1) = e_p(:,:,nxl) |
---|
| 375 | ELSEIF ( inflow_r ) THEN |
---|
| 376 | IF ( .NOT. constant_diffusion ) e_p(:,:,nxr+1) = e_p(:,:,nxr) |
---|
| 377 | ENDIF |
---|
[1] | 378 | |
---|
| 379 | ! |
---|
[1762] | 380 | !-- The same restoration for u at i=nxl and v at j=nys as above must be made |
---|
[1933] | 381 | !-- in case of nest boundaries. This must not be done in case of vertical nesting |
---|
| 382 | !-- mode as in that case the lateral boundaries are actually cyclic. |
---|
| 383 | IF ( nesting_mode /= 'vertical' ) THEN |
---|
| 384 | IF ( nest_bound_s ) THEN |
---|
| 385 | v_p(:,nys,:) = v_p(:,nys-1,:) |
---|
| 386 | ENDIF |
---|
| 387 | IF ( nest_bound_l ) THEN |
---|
| 388 | u_p(:,:,nxl) = u_p(:,:,nxl-1) |
---|
| 389 | ENDIF |
---|
[1762] | 390 | ENDIF |
---|
| 391 | |
---|
| 392 | ! |
---|
[1409] | 393 | !-- Lateral boundary conditions for scalar quantities at the outflow |
---|
| 394 | IF ( outflow_s ) THEN |
---|
| 395 | pt_p(:,nys-1,:) = pt_p(:,nys,:) |
---|
| 396 | IF ( .NOT. constant_diffusion ) e_p(:,nys-1,:) = e_p(:,nys,:) |
---|
[1960] | 397 | IF ( humidity ) THEN |
---|
[1409] | 398 | q_p(:,nys-1,:) = q_p(:,nys,:) |
---|
[1822] | 399 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
[1409] | 400 | qr_p(:,nys-1,:) = qr_p(:,nys,:) |
---|
| 401 | nr_p(:,nys-1,:) = nr_p(:,nys,:) |
---|
[1053] | 402 | ENDIF |
---|
[1409] | 403 | ENDIF |
---|
[1960] | 404 | IF ( passive_scalar ) s_p(:,nys-1,:) = s_p(:,nys,:) |
---|
[1409] | 405 | ELSEIF ( outflow_n ) THEN |
---|
| 406 | pt_p(:,nyn+1,:) = pt_p(:,nyn,:) |
---|
| 407 | IF ( .NOT. constant_diffusion ) e_p(:,nyn+1,:) = e_p(:,nyn,:) |
---|
[1960] | 408 | IF ( humidity ) THEN |
---|
[1409] | 409 | q_p(:,nyn+1,:) = q_p(:,nyn,:) |
---|
[1822] | 410 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
[1409] | 411 | qr_p(:,nyn+1,:) = qr_p(:,nyn,:) |
---|
| 412 | nr_p(:,nyn+1,:) = nr_p(:,nyn,:) |
---|
[1053] | 413 | ENDIF |
---|
[1409] | 414 | ENDIF |
---|
[1960] | 415 | IF ( passive_scalar ) s_p(:,nyn+1,:) = s_p(:,nyn,:) |
---|
[1409] | 416 | ELSEIF ( outflow_l ) THEN |
---|
| 417 | pt_p(:,:,nxl-1) = pt_p(:,:,nxl) |
---|
| 418 | IF ( .NOT. constant_diffusion ) e_p(:,:,nxl-1) = e_p(:,:,nxl) |
---|
[1960] | 419 | IF ( humidity ) THEN |
---|
[1409] | 420 | q_p(:,:,nxl-1) = q_p(:,:,nxl) |
---|
[1822] | 421 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
[1409] | 422 | qr_p(:,:,nxl-1) = qr_p(:,:,nxl) |
---|
| 423 | nr_p(:,:,nxl-1) = nr_p(:,:,nxl) |
---|
[1053] | 424 | ENDIF |
---|
[1409] | 425 | ENDIF |
---|
[1960] | 426 | IF ( passive_scalar ) s_p(:,:,nxl-1) = s_p(:,:,nxl) |
---|
[1409] | 427 | ELSEIF ( outflow_r ) THEN |
---|
| 428 | pt_p(:,:,nxr+1) = pt_p(:,:,nxr) |
---|
| 429 | IF ( .NOT. constant_diffusion ) e_p(:,:,nxr+1) = e_p(:,:,nxr) |
---|
[1960] | 430 | IF ( humidity ) THEN |
---|
[1409] | 431 | q_p(:,:,nxr+1) = q_p(:,:,nxr) |
---|
[1822] | 432 | IF ( cloud_physics .AND. microphysics_seifert ) THEN |
---|
[1409] | 433 | qr_p(:,:,nxr+1) = qr_p(:,:,nxr) |
---|
| 434 | nr_p(:,:,nxr+1) = nr_p(:,:,nxr) |
---|
[1053] | 435 | ENDIF |
---|
[1] | 436 | ENDIF |
---|
[1960] | 437 | IF ( passive_scalar ) s_p(:,:,nxr+1) = s_p(:,:,nxr) |
---|
[1] | 438 | ENDIF |
---|
| 439 | |
---|
| 440 | ! |
---|
[1159] | 441 | !-- Radiation boundary conditions for the velocities at the respective outflow. |
---|
| 442 | !-- The phase velocity is either assumed to the maximum phase velocity that |
---|
| 443 | !-- ensures numerical stability (CFL-condition) or calculated after |
---|
| 444 | !-- Orlanski(1976) and averaged along the outflow boundary. |
---|
[106] | 445 | IF ( outflow_s ) THEN |
---|
[75] | 446 | |
---|
[1159] | 447 | IF ( use_cmax ) THEN |
---|
| 448 | u_p(:,-1,:) = u(:,0,:) |
---|
| 449 | v_p(:,0,:) = v(:,1,:) |
---|
| 450 | w_p(:,-1,:) = w(:,0,:) |
---|
| 451 | ELSEIF ( .NOT. use_cmax ) THEN |
---|
[75] | 452 | |
---|
[978] | 453 | c_max = dy / dt_3d |
---|
[75] | 454 | |
---|
[1353] | 455 | c_u_m_l = 0.0_wp |
---|
| 456 | c_v_m_l = 0.0_wp |
---|
| 457 | c_w_m_l = 0.0_wp |
---|
[978] | 458 | |
---|
[1353] | 459 | c_u_m = 0.0_wp |
---|
| 460 | c_v_m = 0.0_wp |
---|
| 461 | c_w_m = 0.0_wp |
---|
[978] | 462 | |
---|
[75] | 463 | ! |
---|
[996] | 464 | !-- Calculate the phase speeds for u, v, and w, first local and then |
---|
| 465 | !-- average along the outflow boundary. |
---|
| 466 | DO k = nzb+1, nzt+1 |
---|
| 467 | DO i = nxl, nxr |
---|
[75] | 468 | |
---|
[106] | 469 | denom = u_m_s(k,0,i) - u_m_s(k,1,i) |
---|
| 470 | |
---|
[1353] | 471 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 472 | c_u(k,i) = -c_max * ( u(k,0,i) - u_m_s(k,0,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 473 | IF ( c_u(k,i) < 0.0_wp ) THEN |
---|
| 474 | c_u(k,i) = 0.0_wp |
---|
[106] | 475 | ELSEIF ( c_u(k,i) > c_max ) THEN |
---|
| 476 | c_u(k,i) = c_max |
---|
| 477 | ENDIF |
---|
| 478 | ELSE |
---|
| 479 | c_u(k,i) = c_max |
---|
[75] | 480 | ENDIF |
---|
| 481 | |
---|
[106] | 482 | denom = v_m_s(k,1,i) - v_m_s(k,2,i) |
---|
| 483 | |
---|
[1353] | 484 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 485 | c_v(k,i) = -c_max * ( v(k,1,i) - v_m_s(k,1,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 486 | IF ( c_v(k,i) < 0.0_wp ) THEN |
---|
| 487 | c_v(k,i) = 0.0_wp |
---|
[106] | 488 | ELSEIF ( c_v(k,i) > c_max ) THEN |
---|
| 489 | c_v(k,i) = c_max |
---|
| 490 | ENDIF |
---|
| 491 | ELSE |
---|
| 492 | c_v(k,i) = c_max |
---|
[75] | 493 | ENDIF |
---|
| 494 | |
---|
[106] | 495 | denom = w_m_s(k,0,i) - w_m_s(k,1,i) |
---|
[75] | 496 | |
---|
[1353] | 497 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 498 | c_w(k,i) = -c_max * ( w(k,0,i) - w_m_s(k,0,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 499 | IF ( c_w(k,i) < 0.0_wp ) THEN |
---|
| 500 | c_w(k,i) = 0.0_wp |
---|
[106] | 501 | ELSEIF ( c_w(k,i) > c_max ) THEN |
---|
| 502 | c_w(k,i) = c_max |
---|
| 503 | ENDIF |
---|
| 504 | ELSE |
---|
| 505 | c_w(k,i) = c_max |
---|
[75] | 506 | ENDIF |
---|
[106] | 507 | |
---|
[978] | 508 | c_u_m_l(k) = c_u_m_l(k) + c_u(k,i) |
---|
| 509 | c_v_m_l(k) = c_v_m_l(k) + c_v(k,i) |
---|
| 510 | c_w_m_l(k) = c_w_m_l(k) + c_w(k,i) |
---|
[106] | 511 | |
---|
[978] | 512 | ENDDO |
---|
| 513 | ENDDO |
---|
[75] | 514 | |
---|
[978] | 515 | #if defined( __parallel ) |
---|
| 516 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 517 | CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 518 | MPI_SUM, comm1dx, ierr ) |
---|
| 519 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 520 | CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 521 | MPI_SUM, comm1dx, ierr ) |
---|
| 522 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 523 | CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 524 | MPI_SUM, comm1dx, ierr ) |
---|
| 525 | #else |
---|
| 526 | c_u_m = c_u_m_l |
---|
| 527 | c_v_m = c_v_m_l |
---|
| 528 | c_w_m = c_w_m_l |
---|
| 529 | #endif |
---|
| 530 | |
---|
| 531 | c_u_m = c_u_m / (nx+1) |
---|
| 532 | c_v_m = c_v_m / (nx+1) |
---|
| 533 | c_w_m = c_w_m / (nx+1) |
---|
| 534 | |
---|
[75] | 535 | ! |
---|
[978] | 536 | !-- Save old timelevels for the next timestep |
---|
| 537 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
| 538 | u_m_s(:,:,:) = u(:,0:1,:) |
---|
| 539 | v_m_s(:,:,:) = v(:,1:2,:) |
---|
| 540 | w_m_s(:,:,:) = w(:,0:1,:) |
---|
| 541 | ENDIF |
---|
| 542 | |
---|
| 543 | ! |
---|
| 544 | !-- Calculate the new velocities |
---|
[996] | 545 | DO k = nzb+1, nzt+1 |
---|
| 546 | DO i = nxlg, nxrg |
---|
[978] | 547 | u_p(k,-1,i) = u(k,-1,i) - dt_3d * tsc(2) * c_u_m(k) * & |
---|
[75] | 548 | ( u(k,-1,i) - u(k,0,i) ) * ddy |
---|
| 549 | |
---|
[978] | 550 | v_p(k,0,i) = v(k,0,i) - dt_3d * tsc(2) * c_v_m(k) * & |
---|
[106] | 551 | ( v(k,0,i) - v(k,1,i) ) * ddy |
---|
[75] | 552 | |
---|
[978] | 553 | w_p(k,-1,i) = w(k,-1,i) - dt_3d * tsc(2) * c_w_m(k) * & |
---|
[75] | 554 | ( w(k,-1,i) - w(k,0,i) ) * ddy |
---|
[978] | 555 | ENDDO |
---|
[75] | 556 | ENDDO |
---|
| 557 | |
---|
| 558 | ! |
---|
[978] | 559 | !-- Bottom boundary at the outflow |
---|
| 560 | IF ( ibc_uv_b == 0 ) THEN |
---|
[1353] | 561 | u_p(nzb,-1,:) = 0.0_wp |
---|
| 562 | v_p(nzb,0,:) = 0.0_wp |
---|
[978] | 563 | ELSE |
---|
| 564 | u_p(nzb,-1,:) = u_p(nzb+1,-1,:) |
---|
| 565 | v_p(nzb,0,:) = v_p(nzb+1,0,:) |
---|
| 566 | ENDIF |
---|
[1353] | 567 | w_p(nzb,-1,:) = 0.0_wp |
---|
[73] | 568 | |
---|
[75] | 569 | ! |
---|
[978] | 570 | !-- Top boundary at the outflow |
---|
| 571 | IF ( ibc_uv_t == 0 ) THEN |
---|
| 572 | u_p(nzt+1,-1,:) = u_init(nzt+1) |
---|
| 573 | v_p(nzt+1,0,:) = v_init(nzt+1) |
---|
| 574 | ELSE |
---|
[1742] | 575 | u_p(nzt+1,-1,:) = u_p(nzt,-1,:) |
---|
| 576 | v_p(nzt+1,0,:) = v_p(nzt,0,:) |
---|
[978] | 577 | ENDIF |
---|
[1353] | 578 | w_p(nzt:nzt+1,-1,:) = 0.0_wp |
---|
[978] | 579 | |
---|
[75] | 580 | ENDIF |
---|
[73] | 581 | |
---|
[75] | 582 | ENDIF |
---|
[73] | 583 | |
---|
[106] | 584 | IF ( outflow_n ) THEN |
---|
[73] | 585 | |
---|
[1159] | 586 | IF ( use_cmax ) THEN |
---|
| 587 | u_p(:,ny+1,:) = u(:,ny,:) |
---|
| 588 | v_p(:,ny+1,:) = v(:,ny,:) |
---|
| 589 | w_p(:,ny+1,:) = w(:,ny,:) |
---|
| 590 | ELSEIF ( .NOT. use_cmax ) THEN |
---|
[75] | 591 | |
---|
[978] | 592 | c_max = dy / dt_3d |
---|
[75] | 593 | |
---|
[1353] | 594 | c_u_m_l = 0.0_wp |
---|
| 595 | c_v_m_l = 0.0_wp |
---|
| 596 | c_w_m_l = 0.0_wp |
---|
[978] | 597 | |
---|
[1353] | 598 | c_u_m = 0.0_wp |
---|
| 599 | c_v_m = 0.0_wp |
---|
| 600 | c_w_m = 0.0_wp |
---|
[978] | 601 | |
---|
[1] | 602 | ! |
---|
[996] | 603 | !-- Calculate the phase speeds for u, v, and w, first local and then |
---|
| 604 | !-- average along the outflow boundary. |
---|
| 605 | DO k = nzb+1, nzt+1 |
---|
| 606 | DO i = nxl, nxr |
---|
[73] | 607 | |
---|
[106] | 608 | denom = u_m_n(k,ny,i) - u_m_n(k,ny-1,i) |
---|
| 609 | |
---|
[1353] | 610 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 611 | c_u(k,i) = -c_max * ( u(k,ny,i) - u_m_n(k,ny,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 612 | IF ( c_u(k,i) < 0.0_wp ) THEN |
---|
| 613 | c_u(k,i) = 0.0_wp |
---|
[106] | 614 | ELSEIF ( c_u(k,i) > c_max ) THEN |
---|
| 615 | c_u(k,i) = c_max |
---|
| 616 | ENDIF |
---|
| 617 | ELSE |
---|
| 618 | c_u(k,i) = c_max |
---|
[73] | 619 | ENDIF |
---|
| 620 | |
---|
[106] | 621 | denom = v_m_n(k,ny,i) - v_m_n(k,ny-1,i) |
---|
[73] | 622 | |
---|
[1353] | 623 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 624 | c_v(k,i) = -c_max * ( v(k,ny,i) - v_m_n(k,ny,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 625 | IF ( c_v(k,i) < 0.0_wp ) THEN |
---|
| 626 | c_v(k,i) = 0.0_wp |
---|
[106] | 627 | ELSEIF ( c_v(k,i) > c_max ) THEN |
---|
| 628 | c_v(k,i) = c_max |
---|
| 629 | ENDIF |
---|
| 630 | ELSE |
---|
| 631 | c_v(k,i) = c_max |
---|
[73] | 632 | ENDIF |
---|
| 633 | |
---|
[106] | 634 | denom = w_m_n(k,ny,i) - w_m_n(k,ny-1,i) |
---|
[73] | 635 | |
---|
[1353] | 636 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 637 | c_w(k,i) = -c_max * ( w(k,ny,i) - w_m_n(k,ny,i) ) / ( denom * tsc(2) ) |
---|
[1353] | 638 | IF ( c_w(k,i) < 0.0_wp ) THEN |
---|
| 639 | c_w(k,i) = 0.0_wp |
---|
[106] | 640 | ELSEIF ( c_w(k,i) > c_max ) THEN |
---|
| 641 | c_w(k,i) = c_max |
---|
| 642 | ENDIF |
---|
| 643 | ELSE |
---|
| 644 | c_w(k,i) = c_max |
---|
[73] | 645 | ENDIF |
---|
[106] | 646 | |
---|
[978] | 647 | c_u_m_l(k) = c_u_m_l(k) + c_u(k,i) |
---|
| 648 | c_v_m_l(k) = c_v_m_l(k) + c_v(k,i) |
---|
| 649 | c_w_m_l(k) = c_w_m_l(k) + c_w(k,i) |
---|
[106] | 650 | |
---|
[978] | 651 | ENDDO |
---|
| 652 | ENDDO |
---|
[73] | 653 | |
---|
[978] | 654 | #if defined( __parallel ) |
---|
| 655 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 656 | CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 657 | MPI_SUM, comm1dx, ierr ) |
---|
| 658 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 659 | CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 660 | MPI_SUM, comm1dx, ierr ) |
---|
| 661 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dx, ierr ) |
---|
| 662 | CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 663 | MPI_SUM, comm1dx, ierr ) |
---|
| 664 | #else |
---|
| 665 | c_u_m = c_u_m_l |
---|
| 666 | c_v_m = c_v_m_l |
---|
| 667 | c_w_m = c_w_m_l |
---|
| 668 | #endif |
---|
| 669 | |
---|
| 670 | c_u_m = c_u_m / (nx+1) |
---|
| 671 | c_v_m = c_v_m / (nx+1) |
---|
| 672 | c_w_m = c_w_m / (nx+1) |
---|
| 673 | |
---|
[73] | 674 | ! |
---|
[978] | 675 | !-- Save old timelevels for the next timestep |
---|
| 676 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
| 677 | u_m_n(:,:,:) = u(:,ny-1:ny,:) |
---|
| 678 | v_m_n(:,:,:) = v(:,ny-1:ny,:) |
---|
| 679 | w_m_n(:,:,:) = w(:,ny-1:ny,:) |
---|
| 680 | ENDIF |
---|
[73] | 681 | |
---|
[978] | 682 | ! |
---|
| 683 | !-- Calculate the new velocities |
---|
[996] | 684 | DO k = nzb+1, nzt+1 |
---|
| 685 | DO i = nxlg, nxrg |
---|
[978] | 686 | u_p(k,ny+1,i) = u(k,ny+1,i) - dt_3d * tsc(2) * c_u_m(k) * & |
---|
| 687 | ( u(k,ny+1,i) - u(k,ny,i) ) * ddy |
---|
[73] | 688 | |
---|
[978] | 689 | v_p(k,ny+1,i) = v(k,ny+1,i) - dt_3d * tsc(2) * c_v_m(k) * & |
---|
| 690 | ( v(k,ny+1,i) - v(k,ny,i) ) * ddy |
---|
[73] | 691 | |
---|
[978] | 692 | w_p(k,ny+1,i) = w(k,ny+1,i) - dt_3d * tsc(2) * c_w_m(k) * & |
---|
| 693 | ( w(k,ny+1,i) - w(k,ny,i) ) * ddy |
---|
| 694 | ENDDO |
---|
[1] | 695 | ENDDO |
---|
| 696 | |
---|
| 697 | ! |
---|
[978] | 698 | !-- Bottom boundary at the outflow |
---|
| 699 | IF ( ibc_uv_b == 0 ) THEN |
---|
[1353] | 700 | u_p(nzb,ny+1,:) = 0.0_wp |
---|
| 701 | v_p(nzb,ny+1,:) = 0.0_wp |
---|
[978] | 702 | ELSE |
---|
| 703 | u_p(nzb,ny+1,:) = u_p(nzb+1,ny+1,:) |
---|
| 704 | v_p(nzb,ny+1,:) = v_p(nzb+1,ny+1,:) |
---|
| 705 | ENDIF |
---|
[1353] | 706 | w_p(nzb,ny+1,:) = 0.0_wp |
---|
[73] | 707 | |
---|
| 708 | ! |
---|
[978] | 709 | !-- Top boundary at the outflow |
---|
| 710 | IF ( ibc_uv_t == 0 ) THEN |
---|
| 711 | u_p(nzt+1,ny+1,:) = u_init(nzt+1) |
---|
| 712 | v_p(nzt+1,ny+1,:) = v_init(nzt+1) |
---|
| 713 | ELSE |
---|
| 714 | u_p(nzt+1,ny+1,:) = u_p(nzt,nyn+1,:) |
---|
| 715 | v_p(nzt+1,ny+1,:) = v_p(nzt,nyn+1,:) |
---|
| 716 | ENDIF |
---|
[1353] | 717 | w_p(nzt:nzt+1,ny+1,:) = 0.0_wp |
---|
[978] | 718 | |
---|
[1] | 719 | ENDIF |
---|
| 720 | |
---|
[75] | 721 | ENDIF |
---|
| 722 | |
---|
[106] | 723 | IF ( outflow_l ) THEN |
---|
[75] | 724 | |
---|
[1159] | 725 | IF ( use_cmax ) THEN |
---|
[1717] | 726 | u_p(:,:,0) = u(:,:,1) |
---|
| 727 | v_p(:,:,-1) = v(:,:,0) |
---|
[1159] | 728 | w_p(:,:,-1) = w(:,:,0) |
---|
| 729 | ELSEIF ( .NOT. use_cmax ) THEN |
---|
[75] | 730 | |
---|
[978] | 731 | c_max = dx / dt_3d |
---|
[75] | 732 | |
---|
[1353] | 733 | c_u_m_l = 0.0_wp |
---|
| 734 | c_v_m_l = 0.0_wp |
---|
| 735 | c_w_m_l = 0.0_wp |
---|
[978] | 736 | |
---|
[1353] | 737 | c_u_m = 0.0_wp |
---|
| 738 | c_v_m = 0.0_wp |
---|
| 739 | c_w_m = 0.0_wp |
---|
[978] | 740 | |
---|
[1] | 741 | ! |
---|
[996] | 742 | !-- Calculate the phase speeds for u, v, and w, first local and then |
---|
| 743 | !-- average along the outflow boundary. |
---|
| 744 | DO k = nzb+1, nzt+1 |
---|
| 745 | DO j = nys, nyn |
---|
[75] | 746 | |
---|
[106] | 747 | denom = u_m_l(k,j,1) - u_m_l(k,j,2) |
---|
| 748 | |
---|
[1353] | 749 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 750 | c_u(k,j) = -c_max * ( u(k,j,1) - u_m_l(k,j,1) ) / ( denom * tsc(2) ) |
---|
[1353] | 751 | IF ( c_u(k,j) < 0.0_wp ) THEN |
---|
| 752 | c_u(k,j) = 0.0_wp |
---|
[107] | 753 | ELSEIF ( c_u(k,j) > c_max ) THEN |
---|
| 754 | c_u(k,j) = c_max |
---|
[106] | 755 | ENDIF |
---|
| 756 | ELSE |
---|
[107] | 757 | c_u(k,j) = c_max |
---|
[75] | 758 | ENDIF |
---|
| 759 | |
---|
[106] | 760 | denom = v_m_l(k,j,0) - v_m_l(k,j,1) |
---|
[75] | 761 | |
---|
[1353] | 762 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 763 | c_v(k,j) = -c_max * ( v(k,j,0) - v_m_l(k,j,0) ) / ( denom * tsc(2) ) |
---|
[1353] | 764 | IF ( c_v(k,j) < 0.0_wp ) THEN |
---|
| 765 | c_v(k,j) = 0.0_wp |
---|
[106] | 766 | ELSEIF ( c_v(k,j) > c_max ) THEN |
---|
| 767 | c_v(k,j) = c_max |
---|
| 768 | ENDIF |
---|
| 769 | ELSE |
---|
| 770 | c_v(k,j) = c_max |
---|
[75] | 771 | ENDIF |
---|
| 772 | |
---|
[106] | 773 | denom = w_m_l(k,j,0) - w_m_l(k,j,1) |
---|
[75] | 774 | |
---|
[1353] | 775 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 776 | c_w(k,j) = -c_max * ( w(k,j,0) - w_m_l(k,j,0) ) / ( denom * tsc(2) ) |
---|
[1353] | 777 | IF ( c_w(k,j) < 0.0_wp ) THEN |
---|
| 778 | c_w(k,j) = 0.0_wp |
---|
[106] | 779 | ELSEIF ( c_w(k,j) > c_max ) THEN |
---|
| 780 | c_w(k,j) = c_max |
---|
| 781 | ENDIF |
---|
| 782 | ELSE |
---|
| 783 | c_w(k,j) = c_max |
---|
[75] | 784 | ENDIF |
---|
[106] | 785 | |
---|
[978] | 786 | c_u_m_l(k) = c_u_m_l(k) + c_u(k,j) |
---|
| 787 | c_v_m_l(k) = c_v_m_l(k) + c_v(k,j) |
---|
| 788 | c_w_m_l(k) = c_w_m_l(k) + c_w(k,j) |
---|
[106] | 789 | |
---|
[978] | 790 | ENDDO |
---|
| 791 | ENDDO |
---|
[75] | 792 | |
---|
[978] | 793 | #if defined( __parallel ) |
---|
| 794 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 795 | CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 796 | MPI_SUM, comm1dy, ierr ) |
---|
| 797 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 798 | CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 799 | MPI_SUM, comm1dy, ierr ) |
---|
| 800 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 801 | CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 802 | MPI_SUM, comm1dy, ierr ) |
---|
| 803 | #else |
---|
| 804 | c_u_m = c_u_m_l |
---|
| 805 | c_v_m = c_v_m_l |
---|
| 806 | c_w_m = c_w_m_l |
---|
| 807 | #endif |
---|
| 808 | |
---|
| 809 | c_u_m = c_u_m / (ny+1) |
---|
| 810 | c_v_m = c_v_m / (ny+1) |
---|
| 811 | c_w_m = c_w_m / (ny+1) |
---|
| 812 | |
---|
[73] | 813 | ! |
---|
[978] | 814 | !-- Save old timelevels for the next timestep |
---|
| 815 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
| 816 | u_m_l(:,:,:) = u(:,:,1:2) |
---|
| 817 | v_m_l(:,:,:) = v(:,:,0:1) |
---|
| 818 | w_m_l(:,:,:) = w(:,:,0:1) |
---|
| 819 | ENDIF |
---|
| 820 | |
---|
| 821 | ! |
---|
| 822 | !-- Calculate the new velocities |
---|
[996] | 823 | DO k = nzb+1, nzt+1 |
---|
[1113] | 824 | DO j = nysg, nyng |
---|
[978] | 825 | u_p(k,j,0) = u(k,j,0) - dt_3d * tsc(2) * c_u_m(k) * & |
---|
[106] | 826 | ( u(k,j,0) - u(k,j,1) ) * ddx |
---|
[75] | 827 | |
---|
[978] | 828 | v_p(k,j,-1) = v(k,j,-1) - dt_3d * tsc(2) * c_v_m(k) * & |
---|
[75] | 829 | ( v(k,j,-1) - v(k,j,0) ) * ddx |
---|
| 830 | |
---|
[978] | 831 | w_p(k,j,-1) = w(k,j,-1) - dt_3d * tsc(2) * c_w_m(k) * & |
---|
[75] | 832 | ( w(k,j,-1) - w(k,j,0) ) * ddx |
---|
[978] | 833 | ENDDO |
---|
[75] | 834 | ENDDO |
---|
| 835 | |
---|
| 836 | ! |
---|
[978] | 837 | !-- Bottom boundary at the outflow |
---|
| 838 | IF ( ibc_uv_b == 0 ) THEN |
---|
[1353] | 839 | u_p(nzb,:,0) = 0.0_wp |
---|
| 840 | v_p(nzb,:,-1) = 0.0_wp |
---|
[978] | 841 | ELSE |
---|
| 842 | u_p(nzb,:,0) = u_p(nzb+1,:,0) |
---|
| 843 | v_p(nzb,:,-1) = v_p(nzb+1,:,-1) |
---|
| 844 | ENDIF |
---|
[1353] | 845 | w_p(nzb,:,-1) = 0.0_wp |
---|
[1] | 846 | |
---|
[75] | 847 | ! |
---|
[978] | 848 | !-- Top boundary at the outflow |
---|
| 849 | IF ( ibc_uv_t == 0 ) THEN |
---|
[1764] | 850 | u_p(nzt+1,:,0) = u_init(nzt+1) |
---|
[978] | 851 | v_p(nzt+1,:,-1) = v_init(nzt+1) |
---|
| 852 | ELSE |
---|
[1764] | 853 | u_p(nzt+1,:,0) = u_p(nzt,:,0) |
---|
[978] | 854 | v_p(nzt+1,:,-1) = v_p(nzt,:,-1) |
---|
| 855 | ENDIF |
---|
[1353] | 856 | w_p(nzt:nzt+1,:,-1) = 0.0_wp |
---|
[978] | 857 | |
---|
[75] | 858 | ENDIF |
---|
[73] | 859 | |
---|
[75] | 860 | ENDIF |
---|
[73] | 861 | |
---|
[106] | 862 | IF ( outflow_r ) THEN |
---|
[73] | 863 | |
---|
[1159] | 864 | IF ( use_cmax ) THEN |
---|
| 865 | u_p(:,:,nx+1) = u(:,:,nx) |
---|
| 866 | v_p(:,:,nx+1) = v(:,:,nx) |
---|
| 867 | w_p(:,:,nx+1) = w(:,:,nx) |
---|
| 868 | ELSEIF ( .NOT. use_cmax ) THEN |
---|
[75] | 869 | |
---|
[978] | 870 | c_max = dx / dt_3d |
---|
[75] | 871 | |
---|
[1353] | 872 | c_u_m_l = 0.0_wp |
---|
| 873 | c_v_m_l = 0.0_wp |
---|
| 874 | c_w_m_l = 0.0_wp |
---|
[978] | 875 | |
---|
[1353] | 876 | c_u_m = 0.0_wp |
---|
| 877 | c_v_m = 0.0_wp |
---|
| 878 | c_w_m = 0.0_wp |
---|
[978] | 879 | |
---|
[1] | 880 | ! |
---|
[996] | 881 | !-- Calculate the phase speeds for u, v, and w, first local and then |
---|
| 882 | !-- average along the outflow boundary. |
---|
| 883 | DO k = nzb+1, nzt+1 |
---|
| 884 | DO j = nys, nyn |
---|
[73] | 885 | |
---|
[106] | 886 | denom = u_m_r(k,j,nx) - u_m_r(k,j,nx-1) |
---|
| 887 | |
---|
[1353] | 888 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 889 | c_u(k,j) = -c_max * ( u(k,j,nx) - u_m_r(k,j,nx) ) / ( denom * tsc(2) ) |
---|
[1353] | 890 | IF ( c_u(k,j) < 0.0_wp ) THEN |
---|
| 891 | c_u(k,j) = 0.0_wp |
---|
[106] | 892 | ELSEIF ( c_u(k,j) > c_max ) THEN |
---|
| 893 | c_u(k,j) = c_max |
---|
| 894 | ENDIF |
---|
| 895 | ELSE |
---|
| 896 | c_u(k,j) = c_max |
---|
[73] | 897 | ENDIF |
---|
| 898 | |
---|
[106] | 899 | denom = v_m_r(k,j,nx) - v_m_r(k,j,nx-1) |
---|
[73] | 900 | |
---|
[1353] | 901 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 902 | c_v(k,j) = -c_max * ( v(k,j,nx) - v_m_r(k,j,nx) ) / ( denom * tsc(2) ) |
---|
[1353] | 903 | IF ( c_v(k,j) < 0.0_wp ) THEN |
---|
| 904 | c_v(k,j) = 0.0_wp |
---|
[106] | 905 | ELSEIF ( c_v(k,j) > c_max ) THEN |
---|
| 906 | c_v(k,j) = c_max |
---|
| 907 | ENDIF |
---|
| 908 | ELSE |
---|
| 909 | c_v(k,j) = c_max |
---|
[73] | 910 | ENDIF |
---|
| 911 | |
---|
[106] | 912 | denom = w_m_r(k,j,nx) - w_m_r(k,j,nx-1) |
---|
[73] | 913 | |
---|
[1353] | 914 | IF ( denom /= 0.0_wp ) THEN |
---|
[996] | 915 | c_w(k,j) = -c_max * ( w(k,j,nx) - w_m_r(k,j,nx) ) / ( denom * tsc(2) ) |
---|
[1353] | 916 | IF ( c_w(k,j) < 0.0_wp ) THEN |
---|
| 917 | c_w(k,j) = 0.0_wp |
---|
[106] | 918 | ELSEIF ( c_w(k,j) > c_max ) THEN |
---|
| 919 | c_w(k,j) = c_max |
---|
| 920 | ENDIF |
---|
| 921 | ELSE |
---|
| 922 | c_w(k,j) = c_max |
---|
[73] | 923 | ENDIF |
---|
[106] | 924 | |
---|
[978] | 925 | c_u_m_l(k) = c_u_m_l(k) + c_u(k,j) |
---|
| 926 | c_v_m_l(k) = c_v_m_l(k) + c_v(k,j) |
---|
| 927 | c_w_m_l(k) = c_w_m_l(k) + c_w(k,j) |
---|
[106] | 928 | |
---|
[978] | 929 | ENDDO |
---|
| 930 | ENDDO |
---|
[73] | 931 | |
---|
[978] | 932 | #if defined( __parallel ) |
---|
| 933 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 934 | CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 935 | MPI_SUM, comm1dy, ierr ) |
---|
| 936 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 937 | CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 938 | MPI_SUM, comm1dy, ierr ) |
---|
| 939 | IF ( collective_wait ) CALL MPI_BARRIER( comm1dy, ierr ) |
---|
| 940 | CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, & |
---|
| 941 | MPI_SUM, comm1dy, ierr ) |
---|
| 942 | #else |
---|
| 943 | c_u_m = c_u_m_l |
---|
| 944 | c_v_m = c_v_m_l |
---|
| 945 | c_w_m = c_w_m_l |
---|
| 946 | #endif |
---|
| 947 | |
---|
| 948 | c_u_m = c_u_m / (ny+1) |
---|
| 949 | c_v_m = c_v_m / (ny+1) |
---|
| 950 | c_w_m = c_w_m / (ny+1) |
---|
| 951 | |
---|
[73] | 952 | ! |
---|
[978] | 953 | !-- Save old timelevels for the next timestep |
---|
| 954 | IF ( intermediate_timestep_count == 1 ) THEN |
---|
| 955 | u_m_r(:,:,:) = u(:,:,nx-1:nx) |
---|
| 956 | v_m_r(:,:,:) = v(:,:,nx-1:nx) |
---|
| 957 | w_m_r(:,:,:) = w(:,:,nx-1:nx) |
---|
| 958 | ENDIF |
---|
[73] | 959 | |
---|
[978] | 960 | ! |
---|
| 961 | !-- Calculate the new velocities |
---|
[996] | 962 | DO k = nzb+1, nzt+1 |
---|
[1113] | 963 | DO j = nysg, nyng |
---|
[978] | 964 | u_p(k,j,nx+1) = u(k,j,nx+1) - dt_3d * tsc(2) * c_u_m(k) * & |
---|
| 965 | ( u(k,j,nx+1) - u(k,j,nx) ) * ddx |
---|
[73] | 966 | |
---|
[978] | 967 | v_p(k,j,nx+1) = v(k,j,nx+1) - dt_3d * tsc(2) * c_v_m(k) * & |
---|
| 968 | ( v(k,j,nx+1) - v(k,j,nx) ) * ddx |
---|
[73] | 969 | |
---|
[978] | 970 | w_p(k,j,nx+1) = w(k,j,nx+1) - dt_3d * tsc(2) * c_w_m(k) * & |
---|
| 971 | ( w(k,j,nx+1) - w(k,j,nx) ) * ddx |
---|
| 972 | ENDDO |
---|
[73] | 973 | ENDDO |
---|
| 974 | |
---|
| 975 | ! |
---|
[978] | 976 | !-- Bottom boundary at the outflow |
---|
| 977 | IF ( ibc_uv_b == 0 ) THEN |
---|
[1353] | 978 | u_p(nzb,:,nx+1) = 0.0_wp |
---|
| 979 | v_p(nzb,:,nx+1) = 0.0_wp |
---|
[978] | 980 | ELSE |
---|
| 981 | u_p(nzb,:,nx+1) = u_p(nzb+1,:,nx+1) |
---|
| 982 | v_p(nzb,:,nx+1) = v_p(nzb+1,:,nx+1) |
---|
| 983 | ENDIF |
---|
[1353] | 984 | w_p(nzb,:,nx+1) = 0.0_wp |
---|
[73] | 985 | |
---|
| 986 | ! |
---|
[978] | 987 | !-- Top boundary at the outflow |
---|
| 988 | IF ( ibc_uv_t == 0 ) THEN |
---|
| 989 | u_p(nzt+1,:,nx+1) = u_init(nzt+1) |
---|
| 990 | v_p(nzt+1,:,nx+1) = v_init(nzt+1) |
---|
| 991 | ELSE |
---|
| 992 | u_p(nzt+1,:,nx+1) = u_p(nzt,:,nx+1) |
---|
| 993 | v_p(nzt+1,:,nx+1) = v_p(nzt,:,nx+1) |
---|
| 994 | ENDIF |
---|
[1742] | 995 | w_p(nzt:nzt+1,:,nx+1) = 0.0_wp |
---|
[978] | 996 | |
---|
[1] | 997 | ENDIF |
---|
| 998 | |
---|
| 999 | ENDIF |
---|
| 1000 | |
---|
| 1001 | END SUBROUTINE boundary_conds |
---|