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

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

Introduction of nested domain system

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