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

Last change on this file since 1159 was 1159, checked in by fricke, 11 years ago

Bugfix: In case of non-cyclic lateral boundary conditions, Neumann boundary conditions for the velocity components at the outflow are in fact radiation boundary conditions using the maximum phase velocity that ensures numerical stability (CFL-condition).
Logical operator use_cmax is now used instead of bc_lr_dirneu/_neudir.

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