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

Last change on this file since 1114 was 1114, checked in by raasch, 11 years ago

last commit documented

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