source: palm/trunk/SOURCE/boundary_conds.f90 @ 1362

Last change on this file since 1362 was 1362, checked in by hoffmann, 10 years ago

last commit documented

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