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

Last change on this file since 990 was 979, checked in by fricke, 12 years ago

last commit documented

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