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

Last change on this file since 1763 was 1763, checked in by hellstea, 8 years ago

last commit documented

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