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

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

merge fricke branch back into trunk

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