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

Last change on this file since 1398 was 1398, checked in by heinze, 10 years ago

adjustments in case of nudging + bugfix for KIND in CMPLX function

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