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

Last change on this file since 875 was 875, checked in by gryschka, 12 years ago

Bugfix in case of dirichlet inflow bc at the right or north boundary

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