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

Last change on this file since 4086 was 4086, checked in by gronemeier, 5 years ago

Bugfix: use constant-flux layer condition for e in all rans modes (boundary_conds.f90)

  • Property svn:keywords set to Id
File size: 48.2 KB
Line 
1!> @file boundary_conds.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2019 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: boundary_conds.f90 4086 2019-07-11 05:55:44Z gronemeier $
27! Bugfix: use constant-flux layer condition for e in all rans modes
28!
29! 3879 2019-04-08 20:25:23Z knoop
30! Bugfix, do not set boundary conditions for potential temperature in neutral
31! runs.
32!
33! 3655 2019-01-07 16:51:22Z knoop
34! OpenACC port for SPEC
35!
36! 3589 2018-11-30 15:09:51Z suehring
37! Move the control parameter "salsa" from salsa_mod to control_parameters
38! (M. Kurppa)
39!
40! 3582 2018-11-29 19:16:36Z suehring
41! variables documented
42!
43! 3467 2018-10-30 19:05:21Z suehring
44! Implementation of a new aerosol module salsa.
45!
46! 3347 2018-10-15 14:21:08Z suehring
47! Bugfix in setting boundary conditions in offline nesting
48!
49! 3341 2018-10-15 10:31:27Z suehring
50! changes concerning modularization of ocean option
51!
52! 3274 2018-09-24 15:42:55Z knoop
53! Modularization of all bulk cloud physics code components
54!
55! 3241 2018-09-12 15:02:00Z raasch
56! unused variables removed
57!
58! 3183 2018-07-27 14:25:55Z suehring
59! Rename some variables concerning LES-LES as well as offline nesting
60!
61! 3182 2018-07-27 13:36:03Z suehring
62! - Use wall function for e_p and diss_p in case of rans_tke_e
63! - move limitation of diss_p from tcm_prognostic
64!
65! 2938 2018-03-27 15:52:42Z suehring
66! Set boundary condition for TKE and TKE dissipation rate in case of nesting
67! and if parent model operates in RANS mode but child model in LES mode.
68! mode
69!
70! 2793 2018-02-07 10:54:33Z suehring
71! Removed preprocessor directive __chem
72!
73! 2718 2018-01-02 08:49:38Z maronga
74! Corrected "Former revisions" section
75!
76! 2696 2017-12-14 17:12:51Z kanani
77! Change in file header (GPL part)
78! Adjust boundary conditions for e and diss in case of TKE-e closure (TG)
79! Implementation of chemistry module (FK)
80!
81! 2569 2017-10-20 11:54:42Z kanani
82! Removed redundant code for ibc_s_b=1 and ibc_q_b=1
83!
84! 2365 2017-08-21 14:59:59Z kanani
85! Vertical grid nesting implemented: exclude setting vertical velocity to zero
86! on fine grid (SadiqHuq)
87!
88! 2320 2017-07-21 12:47:43Z suehring
89! Remove unused control parameter large_scale_forcing from only-list
90!
91! 2292 2017-06-20 09:51:42Z schwenkel
92! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
93! includes two more prognostic equations for cloud drop concentration (nc) 
94! and cloud water content (qc).
95!
96! 2233 2017-05-30 18:08:54Z suehring
97!
98! 2232 2017-05-30 17:47:52Z suehring
99! Set boundary conditions on topography top using flag method.
100!
101! 2118 2017-01-17 16:38:49Z raasch
102! OpenACC directives removed
103!
104! 2000 2016-08-20 18:09:15Z knoop
105! Forced header and separation lines into 80 columns
106!
107! 1992 2016-08-12 15:14:59Z suehring
108! Adjustments for top boundary condition for passive scalar
109!
110! 1960 2016-07-12 16:34:24Z suehring
111! Treat humidity and passive scalar separately
112!
113! 1823 2016-04-07 08:57:52Z hoffmann
114! Initial version of purely vertical nesting introduced.
115!
116! 1822 2016-04-07 07:49:42Z hoffmann
117! icloud_scheme removed. microphyisics_seifert added.
118!
119! 1764 2016-02-28 12:45:19Z raasch
120! index bug for u_p at left outflow removed
121!
122! 1762 2016-02-25 12:31:13Z hellstea
123! Introduction of nested domain feature
124!
125! 1742 2016-01-13 09:50:06Z raasch
126! bugfix for outflow Neumann boundary conditions at bottom and top
127!
128! 1717 2015-11-11 15:09:47Z raasch
129! Bugfix: index error in outflow conditions for left boundary
130!
131! 1682 2015-10-07 23:56:08Z knoop
132! Code annotations made doxygen readable
133!
134! 1410 2014-05-23 12:16:18Z suehring
135! Bugfix: set dirichlet boundary condition for passive_scalar at model domain
136! top
137!
138! 1399 2014-05-07 11:16:25Z heinze
139! Bugfix: set inflow boundary conditions also if no humidity or passive_scalar
140! is used.
141!
142! 1398 2014-05-07 11:15:00Z heinze
143! Dirichlet-condition at the top for u and v changed to u_init and v_init also
144! for large_scale_forcing
145!
146! 1380 2014-04-28 12:40:45Z heinze
147! Adjust Dirichlet-condition at the top for pt in case of nudging
148!
149! 1361 2014-04-16 15:17:48Z hoffmann
150! Bottom and top boundary conditions of rain water content (qr) and
151! rain drop concentration (nr) changed to Dirichlet
152!
153! 1353 2014-04-08 15:21:23Z heinze
154! REAL constants provided with KIND-attribute
155
156! 1320 2014-03-20 08:40:49Z raasch
157! ONLY-attribute added to USE-statements,
158! kind-parameters added to all INTEGER and REAL declaration statements,
159! kinds are defined in new module kinds,
160! revision history before 2012 removed,
161! comment fields (!:) to be used for variable explanations added to
162! all variable declaration statements
163!
164! 1257 2013-11-08 15:18:40Z raasch
165! loop independent clauses added
166!
167! 1241 2013-10-30 11:36:58Z heinze
168! Adjust ug and vg at each timestep in case of large_scale_forcing
169!
170! 1159 2013-05-21 11:58:22Z fricke
171! Bugfix: Neumann boundary conditions for the velocity components at the
172! outflow are in fact radiation boundary conditions using the maximum phase
173! velocity that ensures numerical stability (CFL-condition).
174! Hence, logical operator use_cmax is now used instead of bc_lr_dirneu/_neudir.
175! Bugfix: In case of use_cmax at the outflow, u, v, w are replaced by
176! u_p, v_p, w_p 
177!
178! 1115 2013-03-26 18:16:16Z hoffmann
179! boundary conditions of two-moment cloud scheme are restricted to Neumann-
180! boundary-conditions
181!
182! 1113 2013-03-10 02:48:14Z raasch
183! GPU-porting
184! dummy argument "range" removed
185! Bugfix: wrong index in loops of radiation boundary condition
186!
187! 1053 2012-11-13 17:11:03Z hoffmann
188! boundary conditions for the two new prognostic equations (nr, qr) of the
189! two-moment cloud scheme
190!
191! 1036 2012-10-22 13:43:42Z raasch
192! code put under GPL (PALM 3.9)
193!
194! 996 2012-09-07 10:41:47Z raasch
195! little reformatting
196!
197! 978 2012-08-09 08:28:32Z fricke
198! Neumann boudnary conditions are added at the inflow boundary for the SGS-TKE.
199! Outflow boundary conditions for the velocity components can be set to Neumann
200! conditions or to radiation conditions with a horizontal averaged phase
201! velocity.
202!
203! 875 2012-04-02 15:35:15Z gryschka
204! Bugfix in case of dirichlet inflow bc at the right or north boundary
205!
206! Revision 1.1  1997/09/12 06:21:34  raasch
207! Initial revision
208!
209!
210! Description:
211! ------------
212!> Boundary conditions for the prognostic quantities.
213!> One additional bottom boundary condition is applied for the TKE (=(u*)**2)
214!> in prandtl_fluxes. The cyclic lateral boundary conditions are implicitly
215!> handled in routine exchange_horiz. Pressure boundary conditions are
216!> explicitly set in routines pres, poisfft, poismg and sor.
217!------------------------------------------------------------------------------!
218 SUBROUTINE boundary_conds
219 
220
221    USE arrays_3d,                                                             &
222        ONLY:  c_u, c_u_m, c_u_m_l, c_v, c_v_m, c_v_m_l, c_w, c_w_m, c_w_m_l,  &
223               diss, diss_p, dzu, e_p, nc_p, nr_p, pt, pt_init, pt_p, q,       &
224               q_p, qc_p, qr_p, s, s_p, sa, sa_p, u, u_init, u_m_l, u_m_n,     &
225               u_m_r, u_m_s, u_p, v, v_init, v_m_l, v_m_n, v_m_r, v_m_s, v_p,  &
226               w, w_p, w_m_l, w_m_n, w_m_r, w_m_s
227
228    USE basic_constants_and_equations_mod,                                     &
229        ONLY:  kappa
230
231    USE bulk_cloud_model_mod,                                                  &
232        ONLY: bulk_cloud_model, microphysics_morrison, microphysics_seifert
233
234    USE chemistry_model_mod,                                                   &
235        ONLY:  chem_boundary_conds
236
237    USE control_parameters,                                                    &
238        ONLY:  air_chemistry, bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r,  &
239               bc_dirichlet_s, bc_radiation_l, bc_radiation_n, bc_radiation_r, &
240               bc_radiation_s, bc_pt_t_val, bc_q_t_val, bc_s_t_val,            &
241               child_domain, constant_diffusion, coupling_mode, dt_3d,         &
242               humidity, ibc_pt_b, ibc_pt_t, ibc_q_b, ibc_q_t, ibc_s_b,        &
243               ibc_s_t, ibc_uv_b, ibc_uv_t, intermediate_timestep_count,       &
244               nesting_offline, neutral, nudging, ocean_mode, passive_scalar,  &
245               rans_mode, rans_tke_e, tsc, salsa, use_cmax
246
247    USE grid_variables,                                                        &
248        ONLY:  ddx, ddy, dx, dy
249
250    USE indices,                                                               &
251        ONLY:  nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzb, nzt
252
253    USE kinds
254
255    USE ocean_mod,                                                             &
256        ONLY:  ibc_sa_t
257
258    USE pegrid
259
260    USE pmc_interface,                                                         &
261        ONLY : nesting_mode, rans_mode_parent
262 
263    USE salsa_mod,                                                             &
264        ONLY:  salsa_boundary_conds       
265
266    USE surface_mod,                                                           &
267        ONLY :  bc_h, surf_def_h, surf_def_v, surf_lsm_h, surf_lsm_v,          &
268                surf_usm_h, surf_usm_v
269
270    USE turbulence_closure_mod,                                                &
271        ONLY:  c_0
272
273    IMPLICIT NONE
274
275    INTEGER(iwp) ::  i  !< grid index x direction
276    INTEGER(iwp) ::  j  !< grid index y direction
277    INTEGER(iwp) ::  k  !< grid index z direction
278    INTEGER(iwp) ::  kb !< variable to set respective boundary value, depends on facing.
279    INTEGER(iwp) ::  l  !< running index boundary type, for up- and downward-facing walls
280    INTEGER(iwp) ::  m  !< running index surface elements
281
282    REAL(wp)    ::  c_max !< maximum phase velocity allowed by CFL criterion, used for outflow boundary condition
283    REAL(wp)    ::  denom !< horizontal gradient of velocity component normal to the outflow boundary
284
285
286!
287!-- Bottom boundary
288    IF ( ibc_uv_b == 1 )  THEN
289       u_p(nzb,:,:) = u_p(nzb+1,:,:)
290       v_p(nzb,:,:) = v_p(nzb+1,:,:)
291    ENDIF
292!
293!-- Set zero vertical velocity at topography top (l=0), or bottom (l=1) in case
294!-- of downward-facing surfaces.
295    DO  l = 0, 1
296!
297!--    Set kb, for upward-facing surfaces value at topography top (k-1) is set,
298!--    for downward-facing surfaces at topography bottom (k+1).
299       kb = MERGE( -1, 1, l == 0 )
300       !$OMP PARALLEL DO PRIVATE( i, j, k )
301       !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
302       !$ACC PRESENT(bc_h, w_p)
303       DO  m = 1, bc_h(l)%ns
304          i = bc_h(l)%i(m)           
305          j = bc_h(l)%j(m)
306          k = bc_h(l)%k(m)
307          w_p(k+kb,j,i) = 0.0_wp
308       ENDDO
309    ENDDO
310
311!
312!-- Top boundary. A nested domain ( ibc_uv_t = 3 ) does not require settings.
313    IF ( ibc_uv_t == 0 )  THEN
314        !$ACC KERNELS PRESENT(u_p, v_p, u_init, v_init)
315        u_p(nzt+1,:,:) = u_init(nzt+1)
316        v_p(nzt+1,:,:) = v_init(nzt+1)
317        !$ACC END KERNELS
318    ELSEIF ( ibc_uv_t == 1 )  THEN
319        u_p(nzt+1,:,:) = u_p(nzt,:,:)
320        v_p(nzt+1,:,:) = v_p(nzt,:,:)
321    ENDIF
322
323!
324!-- Vertical nesting: Vertical velocity not zero at the top of the fine grid
325    IF (  .NOT.  child_domain  .AND.  .NOT.  nesting_offline  .AND.            &
326                 TRIM(coupling_mode) /= 'vnested_fine' )  THEN
327       !$ACC KERNELS PRESENT(w_p)
328       w_p(nzt:nzt+1,:,:) = 0.0_wp  !< nzt is not a prognostic level (but cf. pres)
329       !$ACC END KERNELS
330    ENDIF
331
332!
333!-- Temperature at bottom and top boundary.
334!-- In case of coupled runs (ibc_pt_b = 2) the temperature is given by
335!-- the sea surface temperature of the coupled ocean model.
336!-- Dirichlet
337    IF ( .NOT. neutral )  THEN
338       IF ( ibc_pt_b == 0 )  THEN
339          DO  l = 0, 1
340!     
341!--          Set kb, for upward-facing surfaces value at topography top (k-1)
342!--          is set, for downward-facing surfaces at topography bottom (k+1).
343             kb = MERGE( -1, 1, l == 0 )
344             !$OMP PARALLEL DO PRIVATE( i, j, k )
345             DO  m = 1, bc_h(l)%ns
346                i = bc_h(l)%i(m)           
347                j = bc_h(l)%j(m)
348                k = bc_h(l)%k(m)
349                pt_p(k+kb,j,i) = pt(k+kb,j,i)
350             ENDDO
351          ENDDO
352!     
353!--    Neumann, zero-gradient
354       ELSEIF ( ibc_pt_b == 1 )  THEN
355          DO  l = 0, 1
356!     
357!--          Set kb, for upward-facing surfaces value at topography top (k-1)
358!--          is set, for downward-facing surfaces at topography bottom (k+1).
359             kb = MERGE( -1, 1, l == 0 )
360             !$OMP PARALLEL DO PRIVATE( i, j, k )
361             !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
362             !$ACC PRESENT(bc_h, pt_p)
363             DO  m = 1, bc_h(l)%ns
364                i = bc_h(l)%i(m)           
365                j = bc_h(l)%j(m)
366                k = bc_h(l)%k(m)
367                pt_p(k+kb,j,i) = pt_p(k,j,i)
368             ENDDO
369          ENDDO
370       ENDIF
371       
372!     
373!--    Temperature at top boundary
374       IF ( ibc_pt_t == 0 )  THEN
375           pt_p(nzt+1,:,:) = pt(nzt+1,:,:)
376!     
377!--        In case of nudging adjust top boundary to pt which is
378!--        read in from NUDGING-DATA
379           IF ( nudging )  THEN
380              pt_p(nzt+1,:,:) = pt_init(nzt+1)
381           ENDIF
382       ELSEIF ( ibc_pt_t == 1 )  THEN
383           pt_p(nzt+1,:,:) = pt_p(nzt,:,:)
384       ELSEIF ( ibc_pt_t == 2 )  THEN
385           !$ACC KERNELS PRESENT(pt_p, dzu)
386           pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1)
387           !$ACC END KERNELS
388       ENDIF
389    ENDIF
390
391!
392!-- Boundary conditions for TKE.
393!-- Generally Neumann conditions with de/dz=0 are assumed.
394    IF ( .NOT. constant_diffusion )  THEN
395
396       IF ( .NOT. rans_mode )  THEN
397          DO  l = 0, 1
398!
399!--         Set kb, for upward-facing surfaces value at topography top (k-1) is set,
400!--         for downward-facing surfaces at topography bottom (k+1).
401             kb = MERGE( -1, 1, l == 0 )
402             !$OMP PARALLEL DO PRIVATE( i, j, k )
403             !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
404             !$ACC PRESENT(bc_h, e_p)
405             DO  m = 1, bc_h(l)%ns
406                i = bc_h(l)%i(m)           
407                j = bc_h(l)%j(m)
408                k = bc_h(l)%k(m)
409                e_p(k+kb,j,i) = e_p(k,j,i)
410             ENDDO
411          ENDDO
412       ELSE
413!
414!--       Use wall function within constant-flux layer
415!--       Upward-facing surfaces
416!--       Default surfaces
417          DO  m = 1, surf_def_h(0)%ns
418             i = surf_def_h(0)%i(m)
419             j = surf_def_h(0)%j(m)
420             k = surf_def_h(0)%k(m)
421             e_p(k,j,i) = ( surf_def_h(0)%us(m) / c_0 )**2
422          ENDDO
423!
424!--       Natural surfaces
425          DO  m = 1, surf_lsm_h%ns
426             i = surf_lsm_h%i(m)
427             j = surf_lsm_h%j(m)
428             k = surf_lsm_h%k(m)
429             e_p(k,j,i) = ( surf_lsm_h%us(m) / c_0 )**2
430          ENDDO
431!
432!--       Urban surfaces
433          DO  m = 1, surf_usm_h%ns
434             i = surf_usm_h%i(m)
435             j = surf_usm_h%j(m)
436             k = surf_usm_h%k(m)
437             e_p(k,j,i) = ( surf_usm_h%us(m) / c_0 )**2
438          ENDDO
439!
440!--       Vertical surfaces
441          DO  l = 0, 3
442!
443!--          Default surfaces
444             DO  m = 1, surf_def_v(l)%ns
445                i = surf_def_v(l)%i(m)
446                j = surf_def_v(l)%j(m)
447                k = surf_def_v(l)%k(m)
448                e_p(k,j,i) = ( surf_def_v(l)%us(m) / c_0 )**2
449             ENDDO
450!
451!--          Natural surfaces
452             DO  m = 1, surf_lsm_v(l)%ns
453                i = surf_lsm_v(l)%i(m)
454                j = surf_lsm_v(l)%j(m)
455                k = surf_lsm_v(l)%k(m)
456                e_p(k,j,i) = ( surf_lsm_v(l)%us(m) / c_0 )**2
457             ENDDO
458!
459!--          Urban surfaces
460             DO  m = 1, surf_usm_v(l)%ns
461                i = surf_usm_v(l)%i(m)
462                j = surf_usm_v(l)%j(m)
463                k = surf_usm_v(l)%k(m)
464                e_p(k,j,i) = ( surf_usm_v(l)%us(m) / c_0 )**2
465             ENDDO
466          ENDDO
467       ENDIF
468
469       IF ( .NOT. child_domain )  THEN
470          !$ACC KERNELS PRESENT(e_p)
471          e_p(nzt+1,:,:) = e_p(nzt,:,:)
472          !$ACC END KERNELS
473!
474!--    Nesting case: if parent operates in RANS mode and child in LES mode,
475!--    no TKE is transfered. This case, set Neumann conditions at lateral and
476!--    top child boundaries.
477!--    If not ( both either in RANS or in LES mode ), TKE boundary condition
478!--    is treated in the nesting.
479       ELSE
480
481          IF ( rans_mode_parent  .AND.  .NOT. rans_mode )  THEN
482
483             e_p(nzt+1,:,:) = e_p(nzt,:,:)
484             IF ( bc_dirichlet_l )  e_p(:,:,nxl-1) = e_p(:,:,nxl)
485             IF ( bc_dirichlet_r )  e_p(:,:,nxr+1) = e_p(:,:,nxr)
486             IF ( bc_dirichlet_s )  e_p(:,nys-1,:) = e_p(:,nys,:)
487             IF ( bc_dirichlet_n )  e_p(:,nyn+1,:) = e_p(:,nyn,:)
488
489          ENDIF
490       ENDIF
491    ENDIF
492
493!
494!-- Boundary conditions for TKE dissipation rate.
495    IF ( rans_tke_e )  THEN
496!
497!--    Use wall function within constant-flux layer
498!--    Upward-facing surfaces
499!--    Default surfaces
500       DO  m = 1, surf_def_h(0)%ns
501          i = surf_def_h(0)%i(m)
502          j = surf_def_h(0)%j(m)
503          k = surf_def_h(0)%k(m)
504          diss_p(k,j,i) = surf_def_h(0)%us(m)**3          &
505                        / ( kappa * surf_def_h(0)%z_mo(m) )
506       ENDDO
507!
508!--    Natural surfaces
509       DO  m = 1, surf_lsm_h%ns
510          i = surf_lsm_h%i(m)
511          j = surf_lsm_h%j(m)
512          k = surf_lsm_h%k(m)
513          diss_p(k,j,i) = surf_lsm_h%us(m)**3          &
514                        / ( kappa * surf_lsm_h%z_mo(m) )
515       ENDDO
516!
517!--    Urban surfaces
518       DO  m = 1, surf_usm_h%ns
519          i = surf_usm_h%i(m)
520          j = surf_usm_h%j(m)
521          k = surf_usm_h%k(m)
522          diss_p(k,j,i) = surf_usm_h%us(m)**3          &
523                        / ( kappa * surf_usm_h%z_mo(m) )
524       ENDDO
525!
526!--    Vertical surfaces
527       DO  l = 0, 3
528!
529!--       Default surfaces
530          DO  m = 1, surf_def_v(l)%ns
531             i = surf_def_v(l)%i(m)
532             j = surf_def_v(l)%j(m)
533             k = surf_def_v(l)%k(m)
534             diss_p(k,j,i) = surf_def_v(l)%us(m)**3          &
535                           / ( kappa * surf_def_v(l)%z_mo(m) )
536          ENDDO
537!
538!--       Natural surfaces
539          DO  m = 1, surf_lsm_v(l)%ns
540             i = surf_lsm_v(l)%i(m)
541             j = surf_lsm_v(l)%j(m)
542             k = surf_lsm_v(l)%k(m)
543             diss_p(k,j,i) = surf_lsm_v(l)%us(m)**3          &
544                           / ( kappa * surf_lsm_v(l)%z_mo(m) )
545          ENDDO
546!
547!--       Urban surfaces
548          DO  m = 1, surf_usm_v(l)%ns
549             i = surf_usm_v(l)%i(m)
550             j = surf_usm_v(l)%j(m)
551             k = surf_usm_v(l)%k(m)
552             diss_p(k,j,i) = surf_usm_v(l)%us(m)**3          &
553                           / ( kappa * surf_usm_v(l)%z_mo(m) )
554          ENDDO
555       ENDDO
556!
557!--    Limit change of diss to be between -90% and +100%. Also, set an absolute
558!--    minimum value
559       DO  i = nxl, nxr
560          DO  j = nys, nyn
561             DO  k = nzb, nzt+1
562                diss_p(k,j,i) = MAX( MIN( diss_p(k,j,i),          &
563                                          2.0_wp * diss(k,j,i) ), &
564                                     0.1_wp * diss(k,j,i),        &
565                                     0.0001_wp )
566             ENDDO
567          ENDDO
568       ENDDO
569
570       IF ( .NOT. child_domain )  THEN
571          diss_p(nzt+1,:,:) = diss_p(nzt,:,:)
572       ENDIF
573    ENDIF
574
575!
576!-- Boundary conditions for salinity
577    IF ( ocean_mode )  THEN
578!
579!--    Bottom boundary: Neumann condition because salinity flux is always
580!--    given.
581       DO  l = 0, 1
582!
583!--       Set kb, for upward-facing surfaces value at topography top (k-1) is set,
584!--       for downward-facing surfaces at topography bottom (k+1).
585          kb = MERGE( -1, 1, l == 0 )
586          !$OMP PARALLEL DO PRIVATE( i, j, k )
587          DO  m = 1, bc_h(l)%ns
588             i = bc_h(l)%i(m)           
589             j = bc_h(l)%j(m)
590             k = bc_h(l)%k(m)
591             sa_p(k+kb,j,i) = sa_p(k,j,i)
592          ENDDO
593       ENDDO
594!
595!--    Top boundary: Dirichlet or Neumann
596       IF ( ibc_sa_t == 0 )  THEN
597           sa_p(nzt+1,:,:) = sa(nzt+1,:,:)
598       ELSEIF ( ibc_sa_t == 1 )  THEN
599           sa_p(nzt+1,:,:) = sa_p(nzt,:,:)
600       ENDIF
601
602    ENDIF
603
604!
605!-- Boundary conditions for total water content,
606!-- bottom and top boundary (see also temperature)
607    IF ( humidity )  THEN
608!
609!--    Surface conditions for constant_humidity_flux
610!--    Run loop over all non-natural and natural walls. Note, in wall-datatype
611!--    the k coordinate belongs to the atmospheric grid point, therefore, set
612!--    q_p at k-1
613       IF ( ibc_q_b == 0 ) THEN
614
615          DO  l = 0, 1
616!
617!--          Set kb, for upward-facing surfaces value at topography top (k-1) is set,
618!--          for downward-facing surfaces at topography bottom (k+1).
619             kb = MERGE( -1, 1, l == 0 )
620             !$OMP PARALLEL DO PRIVATE( i, j, k )
621             DO  m = 1, bc_h(l)%ns
622                i = bc_h(l)%i(m)           
623                j = bc_h(l)%j(m)
624                k = bc_h(l)%k(m)
625                q_p(k+kb,j,i) = q(k+kb,j,i)
626             ENDDO
627          ENDDO
628         
629       ELSE
630         
631          DO  l = 0, 1
632!
633!--          Set kb, for upward-facing surfaces value at topography top (k-1) is set,
634!--          for downward-facing surfaces at topography bottom (k+1).
635             kb = MERGE( -1, 1, l == 0 )
636             !$OMP PARALLEL DO PRIVATE( i, j, k )
637             DO  m = 1, bc_h(l)%ns
638                i = bc_h(l)%i(m)           
639                j = bc_h(l)%j(m)
640                k = bc_h(l)%k(m)
641                q_p(k+kb,j,i) = q_p(k,j,i)
642             ENDDO
643          ENDDO
644       ENDIF
645!
646!--    Top boundary
647       IF ( ibc_q_t == 0 ) THEN
648          q_p(nzt+1,:,:) = q(nzt+1,:,:)
649       ELSEIF ( ibc_q_t == 1 ) THEN
650          q_p(nzt+1,:,:) = q_p(nzt,:,:) + bc_q_t_val * dzu(nzt+1)
651       ENDIF
652
653       IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
654!             
655!--       Surface conditions cloud water (Dirichlet)
656!--       Run loop over all non-natural and natural walls. Note, in wall-datatype
657!--       the k coordinate belongs to the atmospheric grid point, therefore, set
658!--       qr_p and nr_p at k-1
659          !$OMP PARALLEL DO PRIVATE( i, j, k )
660          DO  m = 1, bc_h(0)%ns
661             i = bc_h(0)%i(m)           
662             j = bc_h(0)%j(m)
663             k = bc_h(0)%k(m)
664             qc_p(k-1,j,i) = 0.0_wp
665             nc_p(k-1,j,i) = 0.0_wp
666          ENDDO
667!
668!--       Top boundary condition for cloud water (Dirichlet)
669          qc_p(nzt+1,:,:) = 0.0_wp
670          nc_p(nzt+1,:,:) = 0.0_wp
671           
672       ENDIF
673
674       IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
675!             
676!--       Surface conditions rain water (Dirichlet)
677!--       Run loop over all non-natural and natural walls. Note, in wall-datatype
678!--       the k coordinate belongs to the atmospheric grid point, therefore, set
679!--       qr_p and nr_p at k-1
680          !$OMP PARALLEL DO PRIVATE( i, j, k )
681          DO  m = 1, bc_h(0)%ns
682             i = bc_h(0)%i(m)           
683             j = bc_h(0)%j(m)
684             k = bc_h(0)%k(m)
685             qr_p(k-1,j,i) = 0.0_wp
686             nr_p(k-1,j,i) = 0.0_wp
687          ENDDO
688!
689!--       Top boundary condition for rain water (Dirichlet)
690          qr_p(nzt+1,:,:) = 0.0_wp
691          nr_p(nzt+1,:,:) = 0.0_wp
692           
693       ENDIF
694    ENDIF
695!
696!-- Boundary conditions for scalar,
697!-- bottom and top boundary (see also temperature)
698    IF ( passive_scalar )  THEN
699!
700!--    Surface conditions for constant_humidity_flux
701!--    Run loop over all non-natural and natural walls. Note, in wall-datatype
702!--    the k coordinate belongs to the atmospheric grid point, therefore, set
703!--    s_p at k-1
704       IF ( ibc_s_b == 0 ) THEN
705         
706          DO  l = 0, 1
707!
708!--          Set kb, for upward-facing surfaces value at topography top (k-1) is set,
709!--          for downward-facing surfaces at topography bottom (k+1).
710             kb = MERGE( -1, 1, l == 0 )
711             !$OMP PARALLEL DO PRIVATE( i, j, k )
712             DO  m = 1, bc_h(l)%ns
713                i = bc_h(l)%i(m)           
714                j = bc_h(l)%j(m)
715                k = bc_h(l)%k(m)
716                s_p(k+kb,j,i) = s(k+kb,j,i)
717             ENDDO
718          ENDDO
719         
720       ELSE
721         
722          DO  l = 0, 1
723!
724!--          Set kb, for upward-facing surfaces value at topography top (k-1) is set,
725!--          for downward-facing surfaces at topography bottom (k+1).
726             kb = MERGE( -1, 1, l == 0 )
727             !$OMP PARALLEL DO PRIVATE( i, j, k )
728             DO  m = 1, bc_h(l)%ns
729                i = bc_h(l)%i(m)           
730                j = bc_h(l)%j(m)
731                k = bc_h(l)%k(m)
732                s_p(k+kb,j,i) = s_p(k,j,i)
733             ENDDO
734          ENDDO
735       ENDIF
736!
737!--    Top boundary condition
738       IF ( ibc_s_t == 0 )  THEN
739          s_p(nzt+1,:,:) = s(nzt+1,:,:)
740       ELSEIF ( ibc_s_t == 1 )  THEN
741          s_p(nzt+1,:,:) = s_p(nzt,:,:)
742       ELSEIF ( ibc_s_t == 2 )  THEN
743          s_p(nzt+1,:,:) = s_p(nzt,:,:) + bc_s_t_val * dzu(nzt+1)
744       ENDIF
745
746    ENDIF   
747!
748!-- Top/bottom boundary conditions for chemical species
749    IF ( air_chemistry )  CALL chem_boundary_conds( 'set_bc_bottomtop' )
750!
751!-- In case of inflow or nest boundary at the south boundary the boundary for v
752!-- is at nys and in case of inflow or nest boundary at the left boundary the
753!-- boundary for u is at nxl. Since in prognostic_equations (cache optimized
754!-- version) these levels are handled as a prognostic level, boundary values
755!-- have to be restored here.
756!-- For the SGS-TKE, Neumann boundary conditions are used at the inflow.
757    IF ( bc_dirichlet_s )  THEN
758       v_p(:,nys,:) = v_p(:,nys-1,:)
759       IF ( .NOT. constant_diffusion ) e_p(:,nys-1,:) = e_p(:,nys,:)
760    ELSEIF ( bc_dirichlet_n )  THEN
761       IF ( .NOT. constant_diffusion ) e_p(:,nyn+1,:) = e_p(:,nyn,:)
762    ELSEIF ( bc_dirichlet_l ) THEN
763       u_p(:,:,nxl) = u_p(:,:,nxl-1)
764       IF ( .NOT. constant_diffusion ) e_p(:,:,nxl-1) = e_p(:,:,nxl)
765    ELSEIF ( bc_dirichlet_r )  THEN
766       IF ( .NOT. constant_diffusion ) e_p(:,:,nxr+1) = e_p(:,:,nxr)
767    ENDIF
768
769!
770!-- The same restoration for u at i=nxl and v at j=nys as above must be made
771!-- in case of nest boundaries. This must not be done in case of vertical nesting
772!-- mode as in that case the lateral boundaries are actually cyclic.
773    IF ( nesting_mode /= 'vertical'  .OR.  nesting_offline )  THEN
774       IF ( bc_dirichlet_s )  THEN
775          v_p(:,nys,:) = v_p(:,nys-1,:)
776       ENDIF
777       IF ( bc_dirichlet_l )  THEN
778          u_p(:,:,nxl) = u_p(:,:,nxl-1)
779       ENDIF
780    ENDIF
781
782!
783!-- Lateral boundary conditions for scalar quantities at the outflow
784    IF ( bc_radiation_s )  THEN
785       pt_p(:,nys-1,:)     = pt_p(:,nys,:)
786       IF ( .NOT. constant_diffusion )  e_p(:,nys-1,:) = e_p(:,nys,:)
787       IF ( rans_tke_e )  diss_p(:,nys-1,:) = diss_p(:,nys,:)
788       IF ( humidity )  THEN
789          q_p(:,nys-1,:) = q_p(:,nys,:)
790          IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
791             qc_p(:,nys-1,:) = qc_p(:,nys,:)
792             nc_p(:,nys-1,:) = nc_p(:,nys,:)
793          ENDIF
794          IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
795             qr_p(:,nys-1,:) = qr_p(:,nys,:)
796             nr_p(:,nys-1,:) = nr_p(:,nys,:)
797          ENDIF
798       ENDIF
799       IF ( passive_scalar )  s_p(:,nys-1,:) = s_p(:,nys,:)
800    ELSEIF ( bc_radiation_n )  THEN
801       pt_p(:,nyn+1,:)     = pt_p(:,nyn,:)
802       IF ( .NOT. constant_diffusion )  e_p(:,nyn+1,:) = e_p(:,nyn,:)
803       IF ( rans_tke_e )  diss_p(:,nyn+1,:) = diss_p(:,nyn,:)
804       IF ( humidity )  THEN
805          q_p(:,nyn+1,:) = q_p(:,nyn,:)
806          IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
807             qc_p(:,nyn+1,:) = qc_p(:,nyn,:)
808             nc_p(:,nyn+1,:) = nc_p(:,nyn,:)
809          ENDIF
810          IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
811             qr_p(:,nyn+1,:) = qr_p(:,nyn,:)
812             nr_p(:,nyn+1,:) = nr_p(:,nyn,:)
813          ENDIF
814       ENDIF
815       IF ( passive_scalar )  s_p(:,nyn+1,:) = s_p(:,nyn,:)
816    ELSEIF ( bc_radiation_l )  THEN
817       pt_p(:,:,nxl-1)     = pt_p(:,:,nxl)
818       IF ( .NOT. constant_diffusion )  e_p(:,:,nxl-1) = e_p(:,:,nxl)
819       IF ( rans_tke_e )  diss_p(:,:,nxl-1) = diss_p(:,:,nxl)
820       IF ( humidity )  THEN
821          q_p(:,:,nxl-1) = q_p(:,:,nxl)
822          IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
823             qc_p(:,:,nxl-1) = qc_p(:,:,nxl)
824             nc_p(:,:,nxl-1) = nc_p(:,:,nxl)
825          ENDIF
826          IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
827             qr_p(:,:,nxl-1) = qr_p(:,:,nxl)
828             nr_p(:,:,nxl-1) = nr_p(:,:,nxl)
829          ENDIF
830       ENDIF
831       IF ( passive_scalar )  s_p(:,:,nxl-1) = s_p(:,:,nxl)
832    ELSEIF ( bc_radiation_r )  THEN
833       pt_p(:,:,nxr+1)     = pt_p(:,:,nxr)
834       IF ( .NOT. constant_diffusion )  e_p(:,:,nxr+1) = e_p(:,:,nxr)
835       IF ( rans_tke_e )  diss_p(:,:,nxr+1) = diss_p(:,:,nxr)
836       IF ( humidity )  THEN
837          q_p(:,:,nxr+1) = q_p(:,:,nxr)
838          IF ( bulk_cloud_model  .AND.  microphysics_morrison )  THEN
839             qc_p(:,:,nxr+1) = qc_p(:,:,nxr)
840             nc_p(:,:,nxr+1) = nc_p(:,:,nxr)
841          ENDIF
842          IF ( bulk_cloud_model  .AND.  microphysics_seifert )  THEN
843             qr_p(:,:,nxr+1) = qr_p(:,:,nxr)
844             nr_p(:,:,nxr+1) = nr_p(:,:,nxr)
845          ENDIF
846       ENDIF
847       IF ( passive_scalar )  s_p(:,:,nxr+1) = s_p(:,:,nxr)
848    ENDIF
849
850!
851!-- Lateral boundary conditions for chemical species
852    IF ( air_chemistry )  CALL chem_boundary_conds( 'set_bc_lateral' )   
853
854!
855!-- Radiation boundary conditions for the velocities at the respective outflow.
856!-- The phase velocity is either assumed to the maximum phase velocity that
857!-- ensures numerical stability (CFL-condition) or calculated after
858!-- Orlanski(1976) and averaged along the outflow boundary.
859    IF ( bc_radiation_s )  THEN
860
861       IF ( use_cmax )  THEN
862          u_p(:,-1,:) = u(:,0,:)
863          v_p(:,0,:)  = v(:,1,:)
864          w_p(:,-1,:) = w(:,0,:)         
865       ELSEIF ( .NOT. use_cmax )  THEN
866
867          c_max = dy / dt_3d
868
869          c_u_m_l = 0.0_wp 
870          c_v_m_l = 0.0_wp
871          c_w_m_l = 0.0_wp
872
873          c_u_m = 0.0_wp 
874          c_v_m = 0.0_wp
875          c_w_m = 0.0_wp
876
877!
878!--       Calculate the phase speeds for u, v, and w, first local and then
879!--       average along the outflow boundary.
880          DO  k = nzb+1, nzt+1
881             DO  i = nxl, nxr
882
883                denom = u_m_s(k,0,i) - u_m_s(k,1,i)
884
885                IF ( denom /= 0.0_wp )  THEN
886                   c_u(k,i) = -c_max * ( u(k,0,i) - u_m_s(k,0,i) ) / ( denom * tsc(2) )
887                   IF ( c_u(k,i) < 0.0_wp )  THEN
888                      c_u(k,i) = 0.0_wp
889                   ELSEIF ( c_u(k,i) > c_max )  THEN
890                      c_u(k,i) = c_max
891                   ENDIF
892                ELSE
893                   c_u(k,i) = c_max
894                ENDIF
895
896                denom = v_m_s(k,1,i) - v_m_s(k,2,i)
897
898                IF ( denom /= 0.0_wp )  THEN
899                   c_v(k,i) = -c_max * ( v(k,1,i) - v_m_s(k,1,i) ) / ( denom * tsc(2) )
900                   IF ( c_v(k,i) < 0.0_wp )  THEN
901                      c_v(k,i) = 0.0_wp
902                   ELSEIF ( c_v(k,i) > c_max )  THEN
903                      c_v(k,i) = c_max
904                   ENDIF
905                ELSE
906                   c_v(k,i) = c_max
907                ENDIF
908
909                denom = w_m_s(k,0,i) - w_m_s(k,1,i)
910
911                IF ( denom /= 0.0_wp )  THEN
912                   c_w(k,i) = -c_max * ( w(k,0,i) - w_m_s(k,0,i) ) / ( denom * tsc(2) )
913                   IF ( c_w(k,i) < 0.0_wp )  THEN
914                      c_w(k,i) = 0.0_wp
915                   ELSEIF ( c_w(k,i) > c_max )  THEN
916                      c_w(k,i) = c_max
917                   ENDIF
918                ELSE
919                   c_w(k,i) = c_max
920                ENDIF
921
922                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
923                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
924                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
925
926             ENDDO
927          ENDDO
928
929#if defined( __parallel )   
930          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
931          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
932                              MPI_SUM, comm1dx, ierr )   
933          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
934          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
935                              MPI_SUM, comm1dx, ierr ) 
936          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
937          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
938                              MPI_SUM, comm1dx, ierr ) 
939#else
940          c_u_m = c_u_m_l
941          c_v_m = c_v_m_l
942          c_w_m = c_w_m_l
943#endif
944
945          c_u_m = c_u_m / (nx+1)
946          c_v_m = c_v_m / (nx+1)
947          c_w_m = c_w_m / (nx+1)
948
949!
950!--       Save old timelevels for the next timestep
951          IF ( intermediate_timestep_count == 1 )  THEN
952             u_m_s(:,:,:) = u(:,0:1,:)
953             v_m_s(:,:,:) = v(:,1:2,:)
954             w_m_s(:,:,:) = w(:,0:1,:)
955          ENDIF
956
957!
958!--       Calculate the new velocities
959          DO  k = nzb+1, nzt+1
960             DO  i = nxlg, nxrg
961                u_p(k,-1,i) = u(k,-1,i) - dt_3d * tsc(2) * c_u_m(k) *          &
962                                       ( u(k,-1,i) - u(k,0,i) ) * ddy
963
964                v_p(k,0,i)  = v(k,0,i)  - dt_3d * tsc(2) * c_v_m(k) *          &
965                                       ( v(k,0,i) - v(k,1,i) ) * ddy
966
967                w_p(k,-1,i) = w(k,-1,i) - dt_3d * tsc(2) * c_w_m(k) *          &
968                                       ( w(k,-1,i) - w(k,0,i) ) * ddy
969             ENDDO
970          ENDDO
971
972!
973!--       Bottom boundary at the outflow
974          IF ( ibc_uv_b == 0 )  THEN
975             u_p(nzb,-1,:) = 0.0_wp 
976             v_p(nzb,0,:)  = 0.0_wp 
977          ELSE                   
978             u_p(nzb,-1,:) =  u_p(nzb+1,-1,:)
979             v_p(nzb,0,:)  =  v_p(nzb+1,0,:)
980          ENDIF
981          w_p(nzb,-1,:) = 0.0_wp
982
983!
984!--       Top boundary at the outflow
985          IF ( ibc_uv_t == 0 )  THEN
986             u_p(nzt+1,-1,:) = u_init(nzt+1)
987             v_p(nzt+1,0,:)  = v_init(nzt+1)
988          ELSE
989             u_p(nzt+1,-1,:) = u_p(nzt,-1,:)
990             v_p(nzt+1,0,:)  = v_p(nzt,0,:)
991          ENDIF
992          w_p(nzt:nzt+1,-1,:) = 0.0_wp
993
994       ENDIF
995
996    ENDIF
997
998    IF ( bc_radiation_n )  THEN
999
1000       IF ( use_cmax )  THEN
1001          u_p(:,ny+1,:) = u(:,ny,:)
1002          v_p(:,ny+1,:) = v(:,ny,:)
1003          w_p(:,ny+1,:) = w(:,ny,:)         
1004       ELSEIF ( .NOT. use_cmax )  THEN
1005
1006          c_max = dy / dt_3d
1007
1008          c_u_m_l = 0.0_wp 
1009          c_v_m_l = 0.0_wp
1010          c_w_m_l = 0.0_wp
1011
1012          c_u_m = 0.0_wp 
1013          c_v_m = 0.0_wp
1014          c_w_m = 0.0_wp
1015
1016!
1017!--       Calculate the phase speeds for u, v, and w, first local and then
1018!--       average along the outflow boundary.
1019          DO  k = nzb+1, nzt+1
1020             DO  i = nxl, nxr
1021
1022                denom = u_m_n(k,ny,i) - u_m_n(k,ny-1,i)
1023
1024                IF ( denom /= 0.0_wp )  THEN
1025                   c_u(k,i) = -c_max * ( u(k,ny,i) - u_m_n(k,ny,i) ) / ( denom * tsc(2) )
1026                   IF ( c_u(k,i) < 0.0_wp )  THEN
1027                      c_u(k,i) = 0.0_wp
1028                   ELSEIF ( c_u(k,i) > c_max )  THEN
1029                      c_u(k,i) = c_max
1030                   ENDIF
1031                ELSE
1032                   c_u(k,i) = c_max
1033                ENDIF
1034
1035                denom = v_m_n(k,ny,i) - v_m_n(k,ny-1,i)
1036
1037                IF ( denom /= 0.0_wp )  THEN
1038                   c_v(k,i) = -c_max * ( v(k,ny,i) - v_m_n(k,ny,i) ) / ( denom * tsc(2) )
1039                   IF ( c_v(k,i) < 0.0_wp )  THEN
1040                      c_v(k,i) = 0.0_wp
1041                   ELSEIF ( c_v(k,i) > c_max )  THEN
1042                      c_v(k,i) = c_max
1043                   ENDIF
1044                ELSE
1045                   c_v(k,i) = c_max
1046                ENDIF
1047
1048                denom = w_m_n(k,ny,i) - w_m_n(k,ny-1,i)
1049
1050                IF ( denom /= 0.0_wp )  THEN
1051                   c_w(k,i) = -c_max * ( w(k,ny,i) - w_m_n(k,ny,i) ) / ( denom * tsc(2) )
1052                   IF ( c_w(k,i) < 0.0_wp )  THEN
1053                      c_w(k,i) = 0.0_wp
1054                   ELSEIF ( c_w(k,i) > c_max )  THEN
1055                      c_w(k,i) = c_max
1056                   ENDIF
1057                ELSE
1058                   c_w(k,i) = c_max
1059                ENDIF
1060
1061                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
1062                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
1063                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
1064
1065             ENDDO
1066          ENDDO
1067
1068#if defined( __parallel )   
1069          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
1070          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
1071                              MPI_SUM, comm1dx, ierr )   
1072          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
1073          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
1074                              MPI_SUM, comm1dx, ierr ) 
1075          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
1076          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
1077                              MPI_SUM, comm1dx, ierr ) 
1078#else
1079          c_u_m = c_u_m_l
1080          c_v_m = c_v_m_l
1081          c_w_m = c_w_m_l
1082#endif
1083
1084          c_u_m = c_u_m / (nx+1)
1085          c_v_m = c_v_m / (nx+1)
1086          c_w_m = c_w_m / (nx+1)
1087
1088!
1089!--       Save old timelevels for the next timestep
1090          IF ( intermediate_timestep_count == 1 )  THEN
1091                u_m_n(:,:,:) = u(:,ny-1:ny,:)
1092                v_m_n(:,:,:) = v(:,ny-1:ny,:)
1093                w_m_n(:,:,:) = w(:,ny-1:ny,:)
1094          ENDIF
1095
1096!
1097!--       Calculate the new velocities
1098          DO  k = nzb+1, nzt+1
1099             DO  i = nxlg, nxrg
1100                u_p(k,ny+1,i) = u(k,ny+1,i) - dt_3d * tsc(2) * c_u_m(k) *      &
1101                                       ( u(k,ny+1,i) - u(k,ny,i) ) * ddy
1102
1103                v_p(k,ny+1,i) = v(k,ny+1,i)  - dt_3d * tsc(2) * c_v_m(k) *     &
1104                                       ( v(k,ny+1,i) - v(k,ny,i) ) * ddy
1105
1106                w_p(k,ny+1,i) = w(k,ny+1,i) - dt_3d * tsc(2) * c_w_m(k) *      &
1107                                       ( w(k,ny+1,i) - w(k,ny,i) ) * ddy
1108             ENDDO
1109          ENDDO
1110
1111!
1112!--       Bottom boundary at the outflow
1113          IF ( ibc_uv_b == 0 )  THEN
1114             u_p(nzb,ny+1,:) = 0.0_wp
1115             v_p(nzb,ny+1,:) = 0.0_wp   
1116          ELSE                   
1117             u_p(nzb,ny+1,:) =  u_p(nzb+1,ny+1,:)
1118             v_p(nzb,ny+1,:) =  v_p(nzb+1,ny+1,:)
1119          ENDIF
1120          w_p(nzb,ny+1,:) = 0.0_wp
1121
1122!
1123!--       Top boundary at the outflow
1124          IF ( ibc_uv_t == 0 )  THEN
1125             u_p(nzt+1,ny+1,:) = u_init(nzt+1)
1126             v_p(nzt+1,ny+1,:) = v_init(nzt+1)
1127          ELSE
1128             u_p(nzt+1,ny+1,:) = u_p(nzt,nyn+1,:)
1129             v_p(nzt+1,ny+1,:) = v_p(nzt,nyn+1,:)
1130          ENDIF
1131          w_p(nzt:nzt+1,ny+1,:) = 0.0_wp
1132
1133       ENDIF
1134
1135    ENDIF
1136
1137    IF ( bc_radiation_l )  THEN
1138
1139       IF ( use_cmax )  THEN
1140          u_p(:,:,0)  = u(:,:,1)
1141          v_p(:,:,-1) = v(:,:,0)
1142          w_p(:,:,-1) = w(:,:,0)         
1143       ELSEIF ( .NOT. use_cmax )  THEN
1144
1145          c_max = dx / dt_3d
1146
1147          c_u_m_l = 0.0_wp 
1148          c_v_m_l = 0.0_wp
1149          c_w_m_l = 0.0_wp
1150
1151          c_u_m = 0.0_wp 
1152          c_v_m = 0.0_wp
1153          c_w_m = 0.0_wp
1154
1155!
1156!--       Calculate the phase speeds for u, v, and w, first local and then
1157!--       average along the outflow boundary.
1158          DO  k = nzb+1, nzt+1
1159             DO  j = nys, nyn
1160
1161                denom = u_m_l(k,j,1) - u_m_l(k,j,2)
1162
1163                IF ( denom /= 0.0_wp )  THEN
1164                   c_u(k,j) = -c_max * ( u(k,j,1) - u_m_l(k,j,1) ) / ( denom * tsc(2) )
1165                   IF ( c_u(k,j) < 0.0_wp )  THEN
1166                      c_u(k,j) = 0.0_wp
1167                   ELSEIF ( c_u(k,j) > c_max )  THEN
1168                      c_u(k,j) = c_max
1169                   ENDIF
1170                ELSE
1171                   c_u(k,j) = c_max
1172                ENDIF
1173
1174                denom = v_m_l(k,j,0) - v_m_l(k,j,1)
1175
1176                IF ( denom /= 0.0_wp )  THEN
1177                   c_v(k,j) = -c_max * ( v(k,j,0) - v_m_l(k,j,0) ) / ( denom * tsc(2) )
1178                   IF ( c_v(k,j) < 0.0_wp )  THEN
1179                      c_v(k,j) = 0.0_wp
1180                   ELSEIF ( c_v(k,j) > c_max )  THEN
1181                      c_v(k,j) = c_max
1182                   ENDIF
1183                ELSE
1184                   c_v(k,j) = c_max
1185                ENDIF
1186
1187                denom = w_m_l(k,j,0) - w_m_l(k,j,1)
1188
1189                IF ( denom /= 0.0_wp )  THEN
1190                   c_w(k,j) = -c_max * ( w(k,j,0) - w_m_l(k,j,0) ) / ( denom * tsc(2) )
1191                   IF ( c_w(k,j) < 0.0_wp )  THEN
1192                      c_w(k,j) = 0.0_wp
1193                   ELSEIF ( c_w(k,j) > c_max )  THEN
1194                      c_w(k,j) = c_max
1195                   ENDIF
1196                ELSE
1197                   c_w(k,j) = c_max
1198                ENDIF
1199
1200                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
1201                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
1202                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
1203
1204             ENDDO
1205          ENDDO
1206
1207#if defined( __parallel )   
1208          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1209          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
1210                              MPI_SUM, comm1dy, ierr )   
1211          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1212          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
1213                              MPI_SUM, comm1dy, ierr ) 
1214          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1215          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
1216                              MPI_SUM, comm1dy, ierr ) 
1217#else
1218          c_u_m = c_u_m_l
1219          c_v_m = c_v_m_l
1220          c_w_m = c_w_m_l
1221#endif
1222
1223          c_u_m = c_u_m / (ny+1)
1224          c_v_m = c_v_m / (ny+1)
1225          c_w_m = c_w_m / (ny+1)
1226
1227!
1228!--       Save old timelevels for the next timestep
1229          IF ( intermediate_timestep_count == 1 )  THEN
1230                u_m_l(:,:,:) = u(:,:,1:2)
1231                v_m_l(:,:,:) = v(:,:,0:1)
1232                w_m_l(:,:,:) = w(:,:,0:1)
1233          ENDIF
1234
1235!
1236!--       Calculate the new velocities
1237          DO  k = nzb+1, nzt+1
1238             DO  j = nysg, nyng
1239                u_p(k,j,0) = u(k,j,0) - dt_3d * tsc(2) * c_u_m(k) *            &
1240                                       ( u(k,j,0) - u(k,j,1) ) * ddx
1241
1242                v_p(k,j,-1) = v(k,j,-1) - dt_3d * tsc(2) * c_v_m(k) *          &
1243                                       ( v(k,j,-1) - v(k,j,0) ) * ddx
1244
1245                w_p(k,j,-1) = w(k,j,-1) - dt_3d * tsc(2) * c_w_m(k) *          &
1246                                       ( w(k,j,-1) - w(k,j,0) ) * ddx
1247             ENDDO
1248          ENDDO
1249
1250!
1251!--       Bottom boundary at the outflow
1252          IF ( ibc_uv_b == 0 )  THEN
1253             u_p(nzb,:,0)  = 0.0_wp 
1254             v_p(nzb,:,-1) = 0.0_wp
1255          ELSE                   
1256             u_p(nzb,:,0)  =  u_p(nzb+1,:,0)
1257             v_p(nzb,:,-1) =  v_p(nzb+1,:,-1)
1258          ENDIF
1259          w_p(nzb,:,-1) = 0.0_wp
1260
1261!
1262!--       Top boundary at the outflow
1263          IF ( ibc_uv_t == 0 )  THEN
1264             u_p(nzt+1,:,0)  = u_init(nzt+1)
1265             v_p(nzt+1,:,-1) = v_init(nzt+1)
1266          ELSE
1267             u_p(nzt+1,:,0)  = u_p(nzt,:,0)
1268             v_p(nzt+1,:,-1) = v_p(nzt,:,-1)
1269          ENDIF
1270          w_p(nzt:nzt+1,:,-1) = 0.0_wp
1271
1272       ENDIF
1273
1274    ENDIF
1275
1276    IF ( bc_radiation_r )  THEN
1277
1278       IF ( use_cmax )  THEN
1279          u_p(:,:,nx+1) = u(:,:,nx)
1280          v_p(:,:,nx+1) = v(:,:,nx)
1281          w_p(:,:,nx+1) = w(:,:,nx)         
1282       ELSEIF ( .NOT. use_cmax )  THEN
1283
1284          c_max = dx / dt_3d
1285
1286          c_u_m_l = 0.0_wp 
1287          c_v_m_l = 0.0_wp
1288          c_w_m_l = 0.0_wp
1289
1290          c_u_m = 0.0_wp 
1291          c_v_m = 0.0_wp
1292          c_w_m = 0.0_wp
1293
1294!
1295!--       Calculate the phase speeds for u, v, and w, first local and then
1296!--       average along the outflow boundary.
1297          DO  k = nzb+1, nzt+1
1298             DO  j = nys, nyn
1299
1300                denom = u_m_r(k,j,nx) - u_m_r(k,j,nx-1)
1301
1302                IF ( denom /= 0.0_wp )  THEN
1303                   c_u(k,j) = -c_max * ( u(k,j,nx) - u_m_r(k,j,nx) ) / ( denom * tsc(2) )
1304                   IF ( c_u(k,j) < 0.0_wp )  THEN
1305                      c_u(k,j) = 0.0_wp
1306                   ELSEIF ( c_u(k,j) > c_max )  THEN
1307                      c_u(k,j) = c_max
1308                   ENDIF
1309                ELSE
1310                   c_u(k,j) = c_max
1311                ENDIF
1312
1313                denom = v_m_r(k,j,nx) - v_m_r(k,j,nx-1)
1314
1315                IF ( denom /= 0.0_wp )  THEN
1316                   c_v(k,j) = -c_max * ( v(k,j,nx) - v_m_r(k,j,nx) ) / ( denom * tsc(2) )
1317                   IF ( c_v(k,j) < 0.0_wp )  THEN
1318                      c_v(k,j) = 0.0_wp
1319                   ELSEIF ( c_v(k,j) > c_max )  THEN
1320                      c_v(k,j) = c_max
1321                   ENDIF
1322                ELSE
1323                   c_v(k,j) = c_max
1324                ENDIF
1325
1326                denom = w_m_r(k,j,nx) - w_m_r(k,j,nx-1)
1327
1328                IF ( denom /= 0.0_wp )  THEN
1329                   c_w(k,j) = -c_max * ( w(k,j,nx) - w_m_r(k,j,nx) ) / ( denom * tsc(2) )
1330                   IF ( c_w(k,j) < 0.0_wp )  THEN
1331                      c_w(k,j) = 0.0_wp
1332                   ELSEIF ( c_w(k,j) > c_max )  THEN
1333                      c_w(k,j) = c_max
1334                   ENDIF
1335                ELSE
1336                   c_w(k,j) = c_max
1337                ENDIF
1338
1339                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
1340                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
1341                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
1342
1343             ENDDO
1344          ENDDO
1345
1346#if defined( __parallel )   
1347          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1348          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
1349                              MPI_SUM, comm1dy, ierr )   
1350          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1351          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
1352                              MPI_SUM, comm1dy, ierr ) 
1353          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
1354          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
1355                              MPI_SUM, comm1dy, ierr ) 
1356#else
1357          c_u_m = c_u_m_l
1358          c_v_m = c_v_m_l
1359          c_w_m = c_w_m_l
1360#endif
1361
1362          c_u_m = c_u_m / (ny+1)
1363          c_v_m = c_v_m / (ny+1)
1364          c_w_m = c_w_m / (ny+1)
1365
1366!
1367!--       Save old timelevels for the next timestep
1368          IF ( intermediate_timestep_count == 1 )  THEN
1369                u_m_r(:,:,:) = u(:,:,nx-1:nx)
1370                v_m_r(:,:,:) = v(:,:,nx-1:nx)
1371                w_m_r(:,:,:) = w(:,:,nx-1:nx)
1372          ENDIF
1373
1374!
1375!--       Calculate the new velocities
1376          DO  k = nzb+1, nzt+1
1377             DO  j = nysg, nyng
1378                u_p(k,j,nx+1) = u(k,j,nx+1) - dt_3d * tsc(2) * c_u_m(k) *      &
1379                                       ( u(k,j,nx+1) - u(k,j,nx) ) * ddx
1380
1381                v_p(k,j,nx+1) = v(k,j,nx+1) - dt_3d * tsc(2) * c_v_m(k) *      &
1382                                       ( v(k,j,nx+1) - v(k,j,nx) ) * ddx
1383
1384                w_p(k,j,nx+1) = w(k,j,nx+1) - dt_3d * tsc(2) * c_w_m(k) *      &
1385                                       ( w(k,j,nx+1) - w(k,j,nx) ) * ddx
1386             ENDDO
1387          ENDDO
1388
1389!
1390!--       Bottom boundary at the outflow
1391          IF ( ibc_uv_b == 0 )  THEN
1392             u_p(nzb,:,nx+1) = 0.0_wp
1393             v_p(nzb,:,nx+1) = 0.0_wp 
1394          ELSE                   
1395             u_p(nzb,:,nx+1) =  u_p(nzb+1,:,nx+1)
1396             v_p(nzb,:,nx+1) =  v_p(nzb+1,:,nx+1)
1397          ENDIF
1398          w_p(nzb,:,nx+1) = 0.0_wp
1399
1400!
1401!--       Top boundary at the outflow
1402          IF ( ibc_uv_t == 0 )  THEN
1403             u_p(nzt+1,:,nx+1) = u_init(nzt+1)
1404             v_p(nzt+1,:,nx+1) = v_init(nzt+1)
1405          ELSE
1406             u_p(nzt+1,:,nx+1) = u_p(nzt,:,nx+1)
1407             v_p(nzt+1,:,nx+1) = v_p(nzt,:,nx+1)
1408          ENDIF
1409          w_p(nzt:nzt+1,:,nx+1) = 0.0_wp
1410
1411       ENDIF
1412
1413    ENDIF
1414
1415    IF ( salsa )  THEN
1416       CALL salsa_boundary_conds
1417    ENDIF
1418
1419 END SUBROUTINE boundary_conds
Note: See TracBrowser for help on using the repository browser.