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

Last change on this file since 1992 was 1992, checked in by suehring, 8 years ago

Prescribing scalar flux at model top; several bugfixes concering data output of scalars and output of flight data

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