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

Last change on this file since 1113 was 1113, checked in by raasch, 10 years ago

GPU porting of boundary conditions and routine pres; index bug removec from radiation boundary condition

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