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

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

routines for nudging and large scale forcing from external file added

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