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

Last change on this file since 996 was 996, checked in by raasch, 12 years ago

parameter use_prior_plot1d_parameters removed; little reformatting

  • Property svn:keywords set to Id
File size: 26.1 KB
Line 
1 SUBROUTINE boundary_conds( range )
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6! little reformatting
7
8! Former revisions:
9! -----------------
10! $Id: boundary_conds.f90 996 2012-09-07 10:41:47Z raasch $
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 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) ) / ( denom * tsc(2) )
270                   IF ( c_u(k,i) < 0.0 )  THEN
271                      c_u(k,i) = 0.0
272                   ELSEIF ( c_u(k,i) > c_max )  THEN
273                      c_u(k,i) = c_max
274                   ENDIF
275                ELSE
276                   c_u(k,i) = c_max
277                ENDIF
278
279                denom = v_m_s(k,1,i) - v_m_s(k,2,i)
280
281                IF ( denom /= 0.0 )  THEN
282                   c_v(k,i) = -c_max * ( v(k,1,i) - v_m_s(k,1,i) ) / ( denom * tsc(2) )
283                   IF ( c_v(k,i) < 0.0 )  THEN
284                      c_v(k,i) = 0.0
285                   ELSEIF ( c_v(k,i) > c_max )  THEN
286                      c_v(k,i) = c_max
287                   ENDIF
288                ELSE
289                   c_v(k,i) = c_max
290                ENDIF
291
292                denom = w_m_s(k,0,i) - w_m_s(k,1,i)
293
294                IF ( denom /= 0.0 )  THEN
295                   c_w(k,i) = -c_max * ( w(k,0,i) - w_m_s(k,0,i) ) / ( denom * tsc(2) )
296                   IF ( c_w(k,i) < 0.0 )  THEN
297                      c_w(k,i) = 0.0
298                   ELSEIF ( c_w(k,i) > c_max )  THEN
299                      c_w(k,i) = c_max
300                   ENDIF
301                ELSE
302                   c_w(k,i) = c_max
303                ENDIF
304
305                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
306                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
307                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
308
309             ENDDO
310          ENDDO
311
312#if defined( __parallel )   
313          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
314          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
315                              MPI_SUM, comm1dx, ierr )   
316          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
317          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_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_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
321                              MPI_SUM, comm1dx, ierr ) 
322#else
323          c_u_m = c_u_m_l
324          c_v_m = c_v_m_l
325          c_w_m = c_w_m_l
326#endif
327
328          c_u_m = c_u_m / (nx+1)
329          c_v_m = c_v_m / (nx+1)
330          c_w_m = c_w_m / (nx+1)
331
332!
333!--       Save old timelevels for the next timestep
334          IF ( intermediate_timestep_count == 1 )  THEN
335             u_m_s(:,:,:) = u(:,0:1,:)
336             v_m_s(:,:,:) = v(:,1:2,:)
337             w_m_s(:,:,:) = w(:,0:1,:)
338          ENDIF
339
340!
341!--       Calculate the new velocities
342          DO  k = nzb+1, nzt+1
343             DO  i = nxlg, nxrg
344                u_p(k,-1,i) = u(k,-1,i) - dt_3d * tsc(2) * c_u_m(k) *          &
345                                       ( u(k,-1,i) - u(k,0,i) ) * ddy
346
347                v_p(k,0,i)  = v(k,0,i)  - dt_3d * tsc(2) * c_v_m(k) *          &
348                                       ( v(k,0,i) - v(k,1,i) ) * ddy
349
350                w_p(k,-1,i) = w(k,-1,i) - dt_3d * tsc(2) * c_w_m(k) *          &
351                                       ( w(k,-1,i) - w(k,0,i) ) * ddy
352             ENDDO
353          ENDDO
354
355!
356!--       Bottom boundary at the outflow
357          IF ( ibc_uv_b == 0 )  THEN
358             u_p(nzb,-1,:) = 0.0 
359             v_p(nzb,0,:)  = 0.0 
360          ELSE                   
361             u_p(nzb,-1,:) =  u_p(nzb+1,-1,:)
362             v_p(nzb,0,:)  =  v_p(nzb+1,0,:)
363          ENDIF
364          w_p(nzb,-1,:) = 0.0
365
366!
367!--       Top boundary at the outflow
368          IF ( ibc_uv_t == 0 )  THEN
369             u_p(nzt+1,-1,:) = u_init(nzt+1)
370             v_p(nzt+1,0,:)  = v_init(nzt+1)
371          ELSE
372             u_p(nzt+1,-1,:) = u(nzt,-1,:)
373             v_p(nzt+1,0,:)  = v(nzt,0,:)
374          ENDIF
375          w_p(nzt:nzt+1,-1,:) = 0.0
376
377       ENDIF
378
379    ENDIF
380
381    IF ( outflow_n )  THEN
382
383       IF ( bc_ns_neudir )  THEN
384          u(:,ny+1,:) = u(:,ny,:)
385          v(:,ny+1,:) = v(:,ny,:)
386          w(:,ny+1,:) = w(:,ny,:)         
387       ELSEIF ( bc_ns_dirrad )  THEN
388
389          c_max = dy / dt_3d
390
391          c_u_m_l = 0.0 
392          c_v_m_l = 0.0
393          c_w_m_l = 0.0
394
395          c_u_m = 0.0 
396          c_v_m = 0.0
397          c_w_m = 0.0
398
399!
400!--       Calculate the phase speeds for u, v, and w, first local and then
401!--       average along the outflow boundary.
402          DO  k = nzb+1, nzt+1
403             DO  i = nxl, nxr
404
405                denom = u_m_n(k,ny,i) - u_m_n(k,ny-1,i)
406
407                IF ( denom /= 0.0 )  THEN
408                   c_u(k,i) = -c_max * ( u(k,ny,i) - u_m_n(k,ny,i) ) / ( denom * tsc(2) )
409                   IF ( c_u(k,i) < 0.0 )  THEN
410                      c_u(k,i) = 0.0
411                   ELSEIF ( c_u(k,i) > c_max )  THEN
412                      c_u(k,i) = c_max
413                   ENDIF
414                ELSE
415                   c_u(k,i) = c_max
416                ENDIF
417
418                denom = v_m_n(k,ny,i) - v_m_n(k,ny-1,i)
419
420                IF ( denom /= 0.0 )  THEN
421                   c_v(k,i) = -c_max * ( v(k,ny,i) - v_m_n(k,ny,i) ) / ( denom * tsc(2) )
422                   IF ( c_v(k,i) < 0.0 )  THEN
423                      c_v(k,i) = 0.0
424                   ELSEIF ( c_v(k,i) > c_max )  THEN
425                      c_v(k,i) = c_max
426                   ENDIF
427                ELSE
428                   c_v(k,i) = c_max
429                ENDIF
430
431                denom = w_m_n(k,ny,i) - w_m_n(k,ny-1,i)
432
433                IF ( denom /= 0.0 )  THEN
434                   c_w(k,i) = -c_max * ( w(k,ny,i) - w_m_n(k,ny,i) ) / ( denom * tsc(2) )
435                   IF ( c_w(k,i) < 0.0 )  THEN
436                      c_w(k,i) = 0.0
437                   ELSEIF ( c_w(k,i) > c_max )  THEN
438                      c_w(k,i) = c_max
439                   ENDIF
440                ELSE
441                   c_w(k,i) = c_max
442                ENDIF
443
444                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
445                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
446                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
447
448             ENDDO
449          ENDDO
450
451#if defined( __parallel )   
452          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
453          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
454                              MPI_SUM, comm1dx, ierr )   
455          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
456          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
457                              MPI_SUM, comm1dx, ierr ) 
458          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
459          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
460                              MPI_SUM, comm1dx, ierr ) 
461#else
462          c_u_m = c_u_m_l
463          c_v_m = c_v_m_l
464          c_w_m = c_w_m_l
465#endif
466
467          c_u_m = c_u_m / (nx+1)
468          c_v_m = c_v_m / (nx+1)
469          c_w_m = c_w_m / (nx+1)
470
471!
472!--       Save old timelevels for the next timestep
473          IF ( intermediate_timestep_count == 1 )  THEN
474                u_m_n(:,:,:) = u(:,ny-1:ny,:)
475                v_m_n(:,:,:) = v(:,ny-1:ny,:)
476                w_m_n(:,:,:) = w(:,ny-1:ny,:)
477          ENDIF
478
479!
480!--       Calculate the new velocities
481          DO  k = nzb+1, nzt+1
482             DO  i = nxlg, nxrg
483                u_p(k,ny+1,i) = u(k,ny+1,i) - dt_3d * tsc(2) * c_u_m(k) *      &
484                                       ( u(k,ny+1,i) - u(k,ny,i) ) * ddy
485
486                v_p(k,ny+1,i) = v(k,ny+1,i)  - dt_3d * tsc(2) * c_v_m(k) *     &
487                                       ( v(k,ny+1,i) - v(k,ny,i) ) * ddy
488
489                w_p(k,ny+1,i) = w(k,ny+1,i) - dt_3d * tsc(2) * c_w_m(k) *      &
490                                       ( w(k,ny+1,i) - w(k,ny,i) ) * ddy
491             ENDDO
492          ENDDO
493
494!
495!--       Bottom boundary at the outflow
496          IF ( ibc_uv_b == 0 )  THEN
497             u_p(nzb,ny+1,:) = 0.0 
498             v_p(nzb,ny+1,:) = 0.0   
499          ELSE                   
500             u_p(nzb,ny+1,:) =  u_p(nzb+1,ny+1,:)
501             v_p(nzb,ny+1,:) =  v_p(nzb+1,ny+1,:)
502          ENDIF
503          w_p(nzb,ny+1,:) = 0.0
504
505!
506!--       Top boundary at the outflow
507          IF ( ibc_uv_t == 0 )  THEN
508             u_p(nzt+1,ny+1,:) = u_init(nzt+1)
509             v_p(nzt+1,ny+1,:) = v_init(nzt+1)
510          ELSE
511             u_p(nzt+1,ny+1,:) = u_p(nzt,nyn+1,:)
512             v_p(nzt+1,ny+1,:) = v_p(nzt,nyn+1,:)
513          ENDIF
514          w_p(nzt:nzt+1,ny+1,:) = 0.0
515
516       ENDIF
517
518    ENDIF
519
520    IF ( outflow_l )  THEN
521
522       IF ( bc_lr_neudir )  THEN
523          u(:,:,-1) = u(:,:,0)
524          v(:,:,0)  = v(:,:,1)
525          w(:,:,-1) = w(:,:,0)         
526       ELSEIF ( bc_ns_dirrad )  THEN
527
528          c_max = dx / dt_3d
529
530          c_u_m_l = 0.0 
531          c_v_m_l = 0.0
532          c_w_m_l = 0.0
533
534          c_u_m = 0.0 
535          c_v_m = 0.0
536          c_w_m = 0.0
537
538!
539!--       Calculate the phase speeds for u, v, and w, first local and then
540!--       average along the outflow boundary.
541          DO  k = nzb+1, nzt+1
542             DO  j = nys, nyn
543
544                denom = u_m_l(k,j,1) - u_m_l(k,j,2)
545
546                IF ( denom /= 0.0 )  THEN
547                   c_u(k,j) = -c_max * ( u(k,j,1) - u_m_l(k,j,1) ) / ( denom * tsc(2) )
548                   IF ( c_u(k,j) < 0.0 )  THEN
549                      c_u(k,j) = 0.0
550                   ELSEIF ( c_u(k,j) > c_max )  THEN
551                      c_u(k,j) = c_max
552                   ENDIF
553                ELSE
554                   c_u(k,j) = c_max
555                ENDIF
556
557                denom = v_m_l(k,j,0) - v_m_l(k,j,1)
558
559                IF ( denom /= 0.0 )  THEN
560                   c_v(k,j) = -c_max * ( v(k,j,0) - v_m_l(k,j,0) ) / ( denom * tsc(2) )
561                   IF ( c_v(k,j) < 0.0 )  THEN
562                      c_v(k,j) = 0.0
563                   ELSEIF ( c_v(k,j) > c_max )  THEN
564                      c_v(k,j) = c_max
565                   ENDIF
566                ELSE
567                   c_v(k,j) = c_max
568                ENDIF
569
570                denom = w_m_l(k,j,0) - w_m_l(k,j,1)
571
572                IF ( denom /= 0.0 )  THEN
573                   c_w(k,j) = -c_max * ( w(k,j,0) - w_m_l(k,j,0) ) / ( denom * tsc(2) )
574                   IF ( c_w(k,j) < 0.0 )  THEN
575                      c_w(k,j) = 0.0
576                   ELSEIF ( c_w(k,j) > c_max )  THEN
577                      c_w(k,j) = c_max
578                   ENDIF
579                ELSE
580                   c_w(k,j) = c_max
581                ENDIF
582
583                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
584                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
585                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
586
587             ENDDO
588          ENDDO
589
590#if defined( __parallel )   
591          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
592          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
593                              MPI_SUM, comm1dy, ierr )   
594          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
595          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
596                              MPI_SUM, comm1dy, ierr ) 
597          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
598          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
599                              MPI_SUM, comm1dy, ierr ) 
600#else
601          c_u_m = c_u_m_l
602          c_v_m = c_v_m_l
603          c_w_m = c_w_m_l
604#endif
605
606          c_u_m = c_u_m / (ny+1)
607          c_v_m = c_v_m / (ny+1)
608          c_w_m = c_w_m / (ny+1)
609
610!
611!--       Save old timelevels for the next timestep
612          IF ( intermediate_timestep_count == 1 )  THEN
613                u_m_l(:,:,:) = u(:,:,1:2)
614                v_m_l(:,:,:) = v(:,:,0:1)
615                w_m_l(:,:,:) = w(:,:,0:1)
616          ENDIF
617
618!
619!--       Calculate the new velocities
620          DO  k = nzb+1, nzt+1
621             DO  i = nxlg, nxrg
622                u_p(k,j,0) = u(k,j,0) - dt_3d * tsc(2) * c_u_m(k) *            &
623                                       ( u(k,j,0) - u(k,j,1) ) * ddx
624
625                v_p(k,j,-1) = v(k,j,-1) - dt_3d * tsc(2) * c_v_m(k) *          &
626                                       ( v(k,j,-1) - v(k,j,0) ) * ddx
627
628                w_p(k,j,-1) = w(k,j,-1) - dt_3d * tsc(2) * c_w_m(k) *          &
629                                       ( w(k,j,-1) - w(k,j,0) ) * ddx
630             ENDDO
631          ENDDO
632
633!
634!--       Bottom boundary at the outflow
635          IF ( ibc_uv_b == 0 )  THEN
636             u_p(nzb,:,0)  = 0.0 
637             v_p(nzb,:,-1) = 0.0
638          ELSE                   
639             u_p(nzb,:,0)  =  u_p(nzb+1,:,0)
640             v_p(nzb,:,-1) =  v_p(nzb+1,:,-1)
641          ENDIF
642          w_p(nzb,:,-1) = 0.0
643
644!
645!--       Top boundary at the outflow
646          IF ( ibc_uv_t == 0 )  THEN
647             u_p(nzt+1,:,-1) = u_init(nzt+1)
648             v_p(nzt+1,:,-1) = v_init(nzt+1)
649          ELSE
650             u_p(nzt+1,:,-1) = u_p(nzt,:,-1)
651             v_p(nzt+1,:,-1) = v_p(nzt,:,-1)
652          ENDIF
653          w_p(nzt:nzt+1,:,-1) = 0.0
654
655       ENDIF
656
657    ENDIF
658
659    IF ( outflow_r )  THEN
660
661       IF ( bc_lr_dirneu )  THEN
662          u(:,:,nx+1) = u(:,:,nx)
663          v(:,:,nx+1) = v(:,:,nx)
664          w(:,:,nx+1) = w(:,:,nx)         
665       ELSEIF ( bc_ns_dirrad )  THEN
666
667          c_max = dx / dt_3d
668
669          c_u_m_l = 0.0 
670          c_v_m_l = 0.0
671          c_w_m_l = 0.0
672
673          c_u_m = 0.0 
674          c_v_m = 0.0
675          c_w_m = 0.0
676
677!
678!--       Calculate the phase speeds for u, v, and w, first local and then
679!--       average along the outflow boundary.
680          DO  k = nzb+1, nzt+1
681             DO  j = nys, nyn
682
683                denom = u_m_r(k,j,nx) - u_m_r(k,j,nx-1)
684
685                IF ( denom /= 0.0 )  THEN
686                   c_u(k,j) = -c_max * ( u(k,j,nx) - u_m_r(k,j,nx) ) / ( denom * tsc(2) )
687                   IF ( c_u(k,j) < 0.0 )  THEN
688                      c_u(k,j) = 0.0
689                   ELSEIF ( c_u(k,j) > c_max )  THEN
690                      c_u(k,j) = c_max
691                   ENDIF
692                ELSE
693                   c_u(k,j) = c_max
694                ENDIF
695
696                denom = v_m_r(k,j,nx) - v_m_r(k,j,nx-1)
697
698                IF ( denom /= 0.0 )  THEN
699                   c_v(k,j) = -c_max * ( v(k,j,nx) - v_m_r(k,j,nx) ) / ( denom * tsc(2) )
700                   IF ( c_v(k,j) < 0.0 )  THEN
701                      c_v(k,j) = 0.0
702                   ELSEIF ( c_v(k,j) > c_max )  THEN
703                      c_v(k,j) = c_max
704                   ENDIF
705                ELSE
706                   c_v(k,j) = c_max
707                ENDIF
708
709                denom = w_m_r(k,j,nx) - w_m_r(k,j,nx-1)
710
711                IF ( denom /= 0.0 )  THEN
712                   c_w(k,j) = -c_max * ( w(k,j,nx) - w_m_r(k,j,nx) ) / ( denom * tsc(2) )
713                   IF ( c_w(k,j) < 0.0 )  THEN
714                      c_w(k,j) = 0.0
715                   ELSEIF ( c_w(k,j) > c_max )  THEN
716                      c_w(k,j) = c_max
717                   ENDIF
718                ELSE
719                   c_w(k,j) = c_max
720                ENDIF
721
722                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
723                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
724                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
725
726             ENDDO
727          ENDDO
728
729#if defined( __parallel )   
730          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
731          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
732                              MPI_SUM, comm1dy, ierr )   
733          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
734          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
735                              MPI_SUM, comm1dy, ierr ) 
736          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
737          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
738                              MPI_SUM, comm1dy, ierr ) 
739#else
740          c_u_m = c_u_m_l
741          c_v_m = c_v_m_l
742          c_w_m = c_w_m_l
743#endif
744
745          c_u_m = c_u_m / (ny+1)
746          c_v_m = c_v_m / (ny+1)
747          c_w_m = c_w_m / (ny+1)
748
749!
750!--       Save old timelevels for the next timestep
751          IF ( intermediate_timestep_count == 1 )  THEN
752                u_m_r(:,:,:) = u(:,:,nx-1:nx)
753                v_m_r(:,:,:) = v(:,:,nx-1:nx)
754                w_m_r(:,:,:) = w(:,:,nx-1:nx)
755          ENDIF
756
757!
758!--       Calculate the new velocities
759          DO  k = nzb+1, nzt+1
760             DO  i = nxlg, nxrg
761                u_p(k,j,nx+1) = u(k,j,nx+1) - dt_3d * tsc(2) * c_u_m(k) *      &
762                                       ( u(k,j,nx+1) - u(k,j,nx) ) * ddx
763
764                v_p(k,j,nx+1) = v(k,j,nx+1) - dt_3d * tsc(2) * c_v_m(k) *      &
765                                       ( v(k,j,nx+1) - v(k,j,nx) ) * ddx
766
767                w_p(k,j,nx+1) = w(k,j,nx+1) - dt_3d * tsc(2) * c_w_m(k) *      &
768                                       ( w(k,j,nx+1) - w(k,j,nx) ) * ddx
769             ENDDO
770          ENDDO
771
772!
773!--       Bottom boundary at the outflow
774          IF ( ibc_uv_b == 0 )  THEN
775             u_p(nzb,:,nx+1) = 0.0
776             v_p(nzb,:,nx+1) = 0.0 
777          ELSE                   
778             u_p(nzb,:,nx+1) =  u_p(nzb+1,:,nx+1)
779             v_p(nzb,:,nx+1) =  v_p(nzb+1,:,nx+1)
780          ENDIF
781          w_p(nzb,:,nx+1) = 0.0
782
783!
784!--       Top boundary at the outflow
785          IF ( ibc_uv_t == 0 )  THEN
786             u_p(nzt+1,:,nx+1) = u_init(nzt+1)
787             v_p(nzt+1,:,nx+1) = v_init(nzt+1)
788          ELSE
789             u_p(nzt+1,:,nx+1) = u_p(nzt,:,nx+1)
790             v_p(nzt+1,:,nx+1) = v_p(nzt,:,nx+1)
791          ENDIF
792          w(nzt:nzt+1,:,nx+1) = 0.0
793
794       ENDIF
795
796    ENDIF
797
798 END SUBROUTINE boundary_conds
Note: See TracBrowser for help on using the repository browser.