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

Last change on this file since 4277 was 4272, checked in by schwenkel, 5 years ago

further modularization of boundary conditions: moving boundary conditions to their modules

  • Property svn:keywords set to Id
File size: 30.4 KB
RevLine 
[1682]1!> @file boundary_conds.f90
[2000]2!------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[2000]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.
[1036]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!
[3655]17! Copyright 1997-2019 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[484]20! Current revisions:
[1]21! -----------------
[1933]22!
[3589]23!
[1321]24! Former revisions:
25! -----------------
26! $Id: boundary_conds.f90 4272 2019-10-23 15:18:57Z schwenkel $
[4272]27! Further modularization of boundary conditions: moved boundary conditions to
28! respective modules
29!
30! 4268 2019-10-17 11:29:38Z schwenkel
[4268]31! Removing bulk cloud variables to respective module
32!
33! 4182 2019-08-22 15:20:23Z scharf
[4182]34! Corrected "Former revisions" section
35!
36! 4102 2019-07-17 16:00:03Z suehring
[4102]37! - Revise setting for boundary conditions at horizontal walls, use the offset
38!   index that belongs to the data structure instead of pre-calculating
39!   the offset index for each facing.
40! - Set boundary conditions for bulk-cloud quantities also at downward facing
41!   surfaces
42!
43! 4087 2019-07-11 11:35:23Z gronemeier
[4087]44! Add comment
45!
46! 4086 2019-07-11 05:55:44Z gronemeier
[4086]47! Bugfix: use constant-flux layer condition for e in all rans modes
48!
49! 3879 2019-04-08 20:25:23Z knoop
[3717]50! Bugfix, do not set boundary conditions for potential temperature in neutral
51! runs.
52!
53! 3655 2019-01-07 16:51:22Z knoop
[3634]54! OpenACC port for SPEC
[1321]55!
[4182]56! Revision 1.1  1997/09/12 06:21:34  raasch
57! Initial revision
58!
59!
[1]60! Description:
61! ------------
[1682]62!> Boundary conditions for the prognostic quantities.
63!> One additional bottom boundary condition is applied for the TKE (=(u*)**2)
64!> in prandtl_fluxes. The cyclic lateral boundary conditions are implicitly
65!> handled in routine exchange_horiz. Pressure boundary conditions are
66!> explicitly set in routines pres, poisfft, poismg and sor.
[1]67!------------------------------------------------------------------------------!
[1682]68 SUBROUTINE boundary_conds
69 
[1]70
[1320]71    USE arrays_3d,                                                             &
72        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,  &
[4268]73               dzu, pt, pt_init, pt_p, q,                                      &
[4272]74               q_p, s, s_p, u, u_init, u_m_l, u_m_n,                 &
[3241]75               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,  &
76               w, w_p, w_m_l, w_m_n, w_m_r, w_m_s
[2696]77
[1320]78    USE control_parameters,                                                    &
[4272]79        ONLY:  bc_dirichlet_l,                                  &
[3182]80               bc_dirichlet_s, bc_radiation_l, bc_radiation_n, bc_radiation_r, &
81               bc_radiation_s, bc_pt_t_val, bc_q_t_val, bc_s_t_val,            &
[4102]82               child_domain, coupling_mode, dt_3d,                             &
[3582]83               humidity, ibc_pt_b, ibc_pt_t, ibc_q_b, ibc_q_t, ibc_s_b,        &
84               ibc_s_t, ibc_uv_b, ibc_uv_t, intermediate_timestep_count,       &
[4272]85               nesting_offline, neutral, nudging, passive_scalar,  &
86               tsc, use_cmax
[1320]87
88    USE grid_variables,                                                        &
89        ONLY:  ddx, ddy, dx, dy
90
91    USE indices,                                                               &
[3294]92        ONLY:  nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzb, nzt
[1320]93
94    USE kinds
95
[1]96    USE pegrid
97
[1933]98    USE pmc_interface,                                                         &
[4102]99        ONLY : nesting_mode
[3467]100 
[2232]101    USE surface_mod,                                                           &
[4102]102        ONLY :  bc_h
[1933]103
[3129]104    USE turbulence_closure_mod,                                                &
[4102]105        ONLY:  tcm_boundary_conds
[3129]106
[1]107    IMPLICIT NONE
108
[2232]109    INTEGER(iwp) ::  i  !< grid index x direction
110    INTEGER(iwp) ::  j  !< grid index y direction
111    INTEGER(iwp) ::  k  !< grid index z direction
112    INTEGER(iwp) ::  l  !< running index boundary type, for up- and downward-facing walls
113    INTEGER(iwp) ::  m  !< running index surface elements
[1]114
[3562]115    REAL(wp)    ::  c_max !< maximum phase velocity allowed by CFL criterion, used for outflow boundary condition
116    REAL(wp)    ::  denom !< horizontal gradient of velocity component normal to the outflow boundary
[1]117
[73]118
[1]119!
[1113]120!-- Bottom boundary
121    IF ( ibc_uv_b == 1 )  THEN
122       u_p(nzb,:,:) = u_p(nzb+1,:,:)
123       v_p(nzb,:,:) = v_p(nzb+1,:,:)
124    ENDIF
[2232]125!
126!-- Set zero vertical velocity at topography top (l=0), or bottom (l=1) in case
127!-- of downward-facing surfaces.
128    DO  l = 0, 1
129       !$OMP PARALLEL DO PRIVATE( i, j, k )
[3634]130       !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
131       !$ACC PRESENT(bc_h, w_p)
[2232]132       DO  m = 1, bc_h(l)%ns
[4268]133          i = bc_h(l)%i(m)
[2232]134          j = bc_h(l)%j(m)
135          k = bc_h(l)%k(m)
[4102]136          w_p(k+bc_h(l)%koff,j,i) = 0.0_wp
[1113]137       ENDDO
138    ENDDO
139
140!
[1762]141!-- Top boundary. A nested domain ( ibc_uv_t = 3 ) does not require settings.
[1113]142    IF ( ibc_uv_t == 0 )  THEN
[3634]143        !$ACC KERNELS PRESENT(u_p, v_p, u_init, v_init)
[1113]144        u_p(nzt+1,:,:) = u_init(nzt+1)
145        v_p(nzt+1,:,:) = v_init(nzt+1)
[3634]146        !$ACC END KERNELS
[1762]147    ELSEIF ( ibc_uv_t == 1 )  THEN
[1113]148        u_p(nzt+1,:,:) = u_p(nzt,:,:)
149        v_p(nzt+1,:,:) = v_p(nzt,:,:)
150    ENDIF
151
[2365]152!
153!-- Vertical nesting: Vertical velocity not zero at the top of the fine grid
[3347]154    IF (  .NOT.  child_domain  .AND.  .NOT.  nesting_offline  .AND.            &
[2365]155                 TRIM(coupling_mode) /= 'vnested_fine' )  THEN
[3634]156       !$ACC KERNELS PRESENT(w_p)
[2365]157       w_p(nzt:nzt+1,:,:) = 0.0_wp  !< nzt is not a prognostic level (but cf. pres)
[3634]158       !$ACC END KERNELS
[1762]159    ENDIF
160
[1113]161!
[2232]162!-- Temperature at bottom and top boundary.
[1113]163!-- In case of coupled runs (ibc_pt_b = 2) the temperature is given by
164!-- the sea surface temperature of the coupled ocean model.
[2232]165!-- Dirichlet
[3717]166    IF ( .NOT. neutral )  THEN
167       IF ( ibc_pt_b == 0 )  THEN
168          DO  l = 0, 1
169             !$OMP PARALLEL DO PRIVATE( i, j, k )
170             DO  m = 1, bc_h(l)%ns
[4268]171                i = bc_h(l)%i(m)
[3717]172                j = bc_h(l)%j(m)
173                k = bc_h(l)%k(m)
[4102]174                pt_p(k+bc_h(l)%koff,j,i) = pt(k+bc_h(l)%koff,j,i)
[3717]175             ENDDO
[1]176          ENDDO
[3717]177!     
178!--    Neumann, zero-gradient
179       ELSEIF ( ibc_pt_b == 1 )  THEN
180          DO  l = 0, 1
181             !$OMP PARALLEL DO PRIVATE( i, j, k )
182             !$ACC PARALLEL LOOP PRIVATE(i, j, k) &
183             !$ACC PRESENT(bc_h, pt_p)
184             DO  m = 1, bc_h(l)%ns
[4272]185                i = bc_h(l)%i(m)
[3717]186                j = bc_h(l)%j(m)
187                k = bc_h(l)%k(m)
[4102]188                pt_p(k+bc_h(l)%koff,j,i) = pt_p(k,j,i)
[3717]189             ENDDO
[1113]190          ENDDO
[3717]191       ENDIF
192       
193!     
194!--    Temperature at top boundary
195       IF ( ibc_pt_t == 0 )  THEN
196           pt_p(nzt+1,:,:) = pt(nzt+1,:,:)
197!     
198!--        In case of nudging adjust top boundary to pt which is
199!--        read in from NUDGING-DATA
200           IF ( nudging )  THEN
201              pt_p(nzt+1,:,:) = pt_init(nzt+1)
202           ENDIF
203       ELSEIF ( ibc_pt_t == 1 )  THEN
204           pt_p(nzt+1,:,:) = pt_p(nzt,:,:)
205       ELSEIF ( ibc_pt_t == 2 )  THEN
206           !$ACC KERNELS PRESENT(pt_p, dzu)
207           pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1)
208           !$ACC END KERNELS
209       ENDIF
[1113]210    ENDIF
[1]211!
[1960]212!-- Boundary conditions for total water content,
[1113]213!-- bottom and top boundary (see also temperature)
[1960]214    IF ( humidity )  THEN
[1113]215!
216!--    Surface conditions for constant_humidity_flux
[2232]217!--    Run loop over all non-natural and natural walls. Note, in wall-datatype
218!--    the k coordinate belongs to the atmospheric grid point, therefore, set
219!--    q_p at k-1
[1113]220       IF ( ibc_q_b == 0 ) THEN
[2232]221
222          DO  l = 0, 1
223             !$OMP PARALLEL DO PRIVATE( i, j, k )
224             DO  m = 1, bc_h(l)%ns
[4268]225                i = bc_h(l)%i(m)
[2232]226                j = bc_h(l)%j(m)
227                k = bc_h(l)%k(m)
[4102]228                q_p(k+bc_h(l)%koff,j,i) = q(k+bc_h(l)%koff,j,i)
[1]229             ENDDO
230          ENDDO
[4268]231
[1113]232       ELSE
[4268]233
[2232]234          DO  l = 0, 1
235             !$OMP PARALLEL DO PRIVATE( i, j, k )
236             DO  m = 1, bc_h(l)%ns
[4268]237                i = bc_h(l)%i(m)
[2232]238                j = bc_h(l)%j(m)
239                k = bc_h(l)%k(m)
[4102]240                q_p(k+bc_h(l)%koff,j,i) = q_p(k,j,i)
[95]241             ENDDO
242          ENDDO
[1113]243       ENDIF
[95]244!
[1113]245!--    Top boundary
[1462]246       IF ( ibc_q_t == 0 ) THEN
247          q_p(nzt+1,:,:) = q(nzt+1,:,:)
248       ELSEIF ( ibc_q_t == 1 ) THEN
[1992]249          q_p(nzt+1,:,:) = q_p(nzt,:,:) + bc_q_t_val * dzu(nzt+1)
[1462]250       ENDIF
[1409]251    ENDIF
[1]252!
[1960]253!-- Boundary conditions for scalar,
254!-- bottom and top boundary (see also temperature)
255    IF ( passive_scalar )  THEN
256!
257!--    Surface conditions for constant_humidity_flux
[2232]258!--    Run loop over all non-natural and natural walls. Note, in wall-datatype
259!--    the k coordinate belongs to the atmospheric grid point, therefore, set
260!--    s_p at k-1
[1960]261       IF ( ibc_s_b == 0 ) THEN
[4268]262
[2232]263          DO  l = 0, 1
264             !$OMP PARALLEL DO PRIVATE( i, j, k )
265             DO  m = 1, bc_h(l)%ns
[4268]266                i = bc_h(l)%i(m)
[2232]267                j = bc_h(l)%j(m)
268                k = bc_h(l)%k(m)
[4102]269                s_p(k+bc_h(l)%koff,j,i) = s(k+bc_h(l)%koff,j,i)
[1960]270             ENDDO
271          ENDDO
[4268]272
[1960]273       ELSE
[4268]274
[2232]275          DO  l = 0, 1
276             !$OMP PARALLEL DO PRIVATE( i, j, k )
277             DO  m = 1, bc_h(l)%ns
[4268]278                i = bc_h(l)%i(m)
[2232]279                j = bc_h(l)%j(m)
280                k = bc_h(l)%k(m)
[4102]281                s_p(k+bc_h(l)%koff,j,i) = s_p(k,j,i)
[1960]282             ENDDO
283          ENDDO
284       ENDIF
285!
[1992]286!--    Top boundary condition
287       IF ( ibc_s_t == 0 )  THEN
[1960]288          s_p(nzt+1,:,:) = s(nzt+1,:,:)
[1992]289       ELSEIF ( ibc_s_t == 1 )  THEN
290          s_p(nzt+1,:,:) = s_p(nzt,:,:)
291       ELSEIF ( ibc_s_t == 2 )  THEN
292          s_p(nzt+1,:,:) = s_p(nzt,:,:) + bc_s_t_val * dzu(nzt+1)
[1960]293       ENDIF
294
[4102]295    ENDIF 
[1960]296!
[4102]297!-- Set boundary conditions for subgrid TKE and dissipation (RANS mode)
298    CALL tcm_boundary_conds
299!
[1762]300!-- In case of inflow or nest boundary at the south boundary the boundary for v
301!-- is at nys and in case of inflow or nest boundary at the left boundary the
302!-- boundary for u is at nxl. Since in prognostic_equations (cache optimized
303!-- version) these levels are handled as a prognostic level, boundary values
304!-- have to be restored here.
[3182]305    IF ( bc_dirichlet_s )  THEN
[1409]306       v_p(:,nys,:) = v_p(:,nys-1,:)
[3182]307    ELSEIF ( bc_dirichlet_l ) THEN
[1409]308       u_p(:,:,nxl) = u_p(:,:,nxl-1)
309    ENDIF
[1]310
311!
[1762]312!-- The same restoration for u at i=nxl and v at j=nys as above must be made
[1933]313!-- in case of nest boundaries. This must not be done in case of vertical nesting
[3182]314!-- mode as in that case the lateral boundaries are actually cyclic.
[4102]315!-- Lateral oundary conditions for TKE and dissipation are set
316!-- in tcm_boundary_conds.
[3182]317    IF ( nesting_mode /= 'vertical'  .OR.  nesting_offline )  THEN
318       IF ( bc_dirichlet_s )  THEN
[1933]319          v_p(:,nys,:) = v_p(:,nys-1,:)
320       ENDIF
[3182]321       IF ( bc_dirichlet_l )  THEN
[1933]322          u_p(:,:,nxl) = u_p(:,:,nxl-1)
323       ENDIF
[1762]324    ENDIF
325
326!
[4102]327!-- Lateral boundary conditions for scalar quantities at the outflow.
328!-- Lateral oundary conditions for TKE and dissipation are set
329!-- in tcm_boundary_conds.
[3182]330    IF ( bc_radiation_s )  THEN
[1409]331       pt_p(:,nys-1,:)     = pt_p(:,nys,:)
[1960]332       IF ( humidity )  THEN
[1409]333          q_p(:,nys-1,:) = q_p(:,nys,:)
334       ENDIF
[1960]335       IF ( passive_scalar )  s_p(:,nys-1,:) = s_p(:,nys,:)
[3182]336    ELSEIF ( bc_radiation_n )  THEN
[1409]337       pt_p(:,nyn+1,:)     = pt_p(:,nyn,:)
[1960]338       IF ( humidity )  THEN
[1409]339          q_p(:,nyn+1,:) = q_p(:,nyn,:)
340       ENDIF
[1960]341       IF ( passive_scalar )  s_p(:,nyn+1,:) = s_p(:,nyn,:)
[3182]342    ELSEIF ( bc_radiation_l )  THEN
[1409]343       pt_p(:,:,nxl-1)     = pt_p(:,:,nxl)
[1960]344       IF ( humidity )  THEN
[1409]345          q_p(:,:,nxl-1) = q_p(:,:,nxl)
346       ENDIF
[1960]347       IF ( passive_scalar )  s_p(:,:,nxl-1) = s_p(:,:,nxl)
[3182]348    ELSEIF ( bc_radiation_r )  THEN
[1409]349       pt_p(:,:,nxr+1)     = pt_p(:,:,nxr)
[1960]350       IF ( humidity )  THEN
[1409]351          q_p(:,:,nxr+1) = q_p(:,:,nxr)
[1]352       ENDIF
[1960]353       IF ( passive_scalar )  s_p(:,:,nxr+1) = s_p(:,:,nxr)
[1]354    ENDIF
355
356!
[1159]357!-- Radiation boundary conditions for the velocities at the respective outflow.
358!-- The phase velocity is either assumed to the maximum phase velocity that
359!-- ensures numerical stability (CFL-condition) or calculated after
360!-- Orlanski(1976) and averaged along the outflow boundary.
[3182]361    IF ( bc_radiation_s )  THEN
[75]362
[1159]363       IF ( use_cmax )  THEN
364          u_p(:,-1,:) = u(:,0,:)
365          v_p(:,0,:)  = v(:,1,:)
[4268]366          w_p(:,-1,:) = w(:,0,:)
[1159]367       ELSEIF ( .NOT. use_cmax )  THEN
[75]368
[978]369          c_max = dy / dt_3d
[75]370
[1353]371          c_u_m_l = 0.0_wp 
372          c_v_m_l = 0.0_wp
373          c_w_m_l = 0.0_wp
[978]374
[1353]375          c_u_m = 0.0_wp 
376          c_v_m = 0.0_wp
377          c_w_m = 0.0_wp
[978]378
[75]379!
[996]380!--       Calculate the phase speeds for u, v, and w, first local and then
381!--       average along the outflow boundary.
382          DO  k = nzb+1, nzt+1
383             DO  i = nxl, nxr
[75]384
[106]385                denom = u_m_s(k,0,i) - u_m_s(k,1,i)
386
[1353]387                IF ( denom /= 0.0_wp )  THEN
[996]388                   c_u(k,i) = -c_max * ( u(k,0,i) - u_m_s(k,0,i) ) / ( denom * tsc(2) )
[1353]389                   IF ( c_u(k,i) < 0.0_wp )  THEN
390                      c_u(k,i) = 0.0_wp
[106]391                   ELSEIF ( c_u(k,i) > c_max )  THEN
392                      c_u(k,i) = c_max
393                   ENDIF
394                ELSE
395                   c_u(k,i) = c_max
[75]396                ENDIF
397
[106]398                denom = v_m_s(k,1,i) - v_m_s(k,2,i)
399
[1353]400                IF ( denom /= 0.0_wp )  THEN
[996]401                   c_v(k,i) = -c_max * ( v(k,1,i) - v_m_s(k,1,i) ) / ( denom * tsc(2) )
[1353]402                   IF ( c_v(k,i) < 0.0_wp )  THEN
403                      c_v(k,i) = 0.0_wp
[106]404                   ELSEIF ( c_v(k,i) > c_max )  THEN
405                      c_v(k,i) = c_max
406                   ENDIF
407                ELSE
408                   c_v(k,i) = c_max
[75]409                ENDIF
410
[106]411                denom = w_m_s(k,0,i) - w_m_s(k,1,i)
[75]412
[1353]413                IF ( denom /= 0.0_wp )  THEN
[996]414                   c_w(k,i) = -c_max * ( w(k,0,i) - w_m_s(k,0,i) ) / ( denom * tsc(2) )
[1353]415                   IF ( c_w(k,i) < 0.0_wp )  THEN
416                      c_w(k,i) = 0.0_wp
[106]417                   ELSEIF ( c_w(k,i) > c_max )  THEN
418                      c_w(k,i) = c_max
419                   ENDIF
420                ELSE
421                   c_w(k,i) = c_max
[75]422                ENDIF
[106]423
[978]424                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
425                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
426                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
[106]427
[978]428             ENDDO
429          ENDDO
[75]430
[4268]431#if defined( __parallel )
[978]432          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
433          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
434                              MPI_SUM, comm1dx, ierr )   
435          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
436          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
437                              MPI_SUM, comm1dx, ierr ) 
438          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
439          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
440                              MPI_SUM, comm1dx, ierr ) 
441#else
442          c_u_m = c_u_m_l
443          c_v_m = c_v_m_l
444          c_w_m = c_w_m_l
445#endif
446
447          c_u_m = c_u_m / (nx+1)
448          c_v_m = c_v_m / (nx+1)
449          c_w_m = c_w_m / (nx+1)
450
[75]451!
[978]452!--       Save old timelevels for the next timestep
453          IF ( intermediate_timestep_count == 1 )  THEN
454             u_m_s(:,:,:) = u(:,0:1,:)
455             v_m_s(:,:,:) = v(:,1:2,:)
456             w_m_s(:,:,:) = w(:,0:1,:)
457          ENDIF
458
459!
460!--       Calculate the new velocities
[996]461          DO  k = nzb+1, nzt+1
462             DO  i = nxlg, nxrg
[978]463                u_p(k,-1,i) = u(k,-1,i) - dt_3d * tsc(2) * c_u_m(k) *          &
[75]464                                       ( u(k,-1,i) - u(k,0,i) ) * ddy
465
[978]466                v_p(k,0,i)  = v(k,0,i)  - dt_3d * tsc(2) * c_v_m(k) *          &
[106]467                                       ( v(k,0,i) - v(k,1,i) ) * ddy
[75]468
[978]469                w_p(k,-1,i) = w(k,-1,i) - dt_3d * tsc(2) * c_w_m(k) *          &
[75]470                                       ( w(k,-1,i) - w(k,0,i) ) * ddy
[978]471             ENDDO
[75]472          ENDDO
473
474!
[978]475!--       Bottom boundary at the outflow
476          IF ( ibc_uv_b == 0 )  THEN
[1353]477             u_p(nzb,-1,:) = 0.0_wp 
478             v_p(nzb,0,:)  = 0.0_wp 
[4268]479          ELSE
[978]480             u_p(nzb,-1,:) =  u_p(nzb+1,-1,:)
481             v_p(nzb,0,:)  =  v_p(nzb+1,0,:)
482          ENDIF
[1353]483          w_p(nzb,-1,:) = 0.0_wp
[73]484
[75]485!
[978]486!--       Top boundary at the outflow
487          IF ( ibc_uv_t == 0 )  THEN
488             u_p(nzt+1,-1,:) = u_init(nzt+1)
489             v_p(nzt+1,0,:)  = v_init(nzt+1)
490          ELSE
[1742]491             u_p(nzt+1,-1,:) = u_p(nzt,-1,:)
492             v_p(nzt+1,0,:)  = v_p(nzt,0,:)
[978]493          ENDIF
[1353]494          w_p(nzt:nzt+1,-1,:) = 0.0_wp
[978]495
[75]496       ENDIF
[73]497
[75]498    ENDIF
[73]499
[3182]500    IF ( bc_radiation_n )  THEN
[73]501
[1159]502       IF ( use_cmax )  THEN
503          u_p(:,ny+1,:) = u(:,ny,:)
504          v_p(:,ny+1,:) = v(:,ny,:)
[4268]505          w_p(:,ny+1,:) = w(:,ny,:)
[1159]506       ELSEIF ( .NOT. use_cmax )  THEN
[75]507
[978]508          c_max = dy / dt_3d
[75]509
[1353]510          c_u_m_l = 0.0_wp 
511          c_v_m_l = 0.0_wp
512          c_w_m_l = 0.0_wp
[978]513
[1353]514          c_u_m = 0.0_wp 
515          c_v_m = 0.0_wp
516          c_w_m = 0.0_wp
[978]517
[1]518!
[996]519!--       Calculate the phase speeds for u, v, and w, first local and then
520!--       average along the outflow boundary.
521          DO  k = nzb+1, nzt+1
522             DO  i = nxl, nxr
[73]523
[106]524                denom = u_m_n(k,ny,i) - u_m_n(k,ny-1,i)
525
[1353]526                IF ( denom /= 0.0_wp )  THEN
[996]527                   c_u(k,i) = -c_max * ( u(k,ny,i) - u_m_n(k,ny,i) ) / ( denom * tsc(2) )
[1353]528                   IF ( c_u(k,i) < 0.0_wp )  THEN
529                      c_u(k,i) = 0.0_wp
[106]530                   ELSEIF ( c_u(k,i) > c_max )  THEN
531                      c_u(k,i) = c_max
532                   ENDIF
533                ELSE
534                   c_u(k,i) = c_max
[73]535                ENDIF
536
[106]537                denom = v_m_n(k,ny,i) - v_m_n(k,ny-1,i)
[73]538
[1353]539                IF ( denom /= 0.0_wp )  THEN
[996]540                   c_v(k,i) = -c_max * ( v(k,ny,i) - v_m_n(k,ny,i) ) / ( denom * tsc(2) )
[1353]541                   IF ( c_v(k,i) < 0.0_wp )  THEN
542                      c_v(k,i) = 0.0_wp
[106]543                   ELSEIF ( c_v(k,i) > c_max )  THEN
544                      c_v(k,i) = c_max
545                   ENDIF
546                ELSE
547                   c_v(k,i) = c_max
[73]548                ENDIF
549
[106]550                denom = w_m_n(k,ny,i) - w_m_n(k,ny-1,i)
[73]551
[1353]552                IF ( denom /= 0.0_wp )  THEN
[996]553                   c_w(k,i) = -c_max * ( w(k,ny,i) - w_m_n(k,ny,i) ) / ( denom * tsc(2) )
[1353]554                   IF ( c_w(k,i) < 0.0_wp )  THEN
555                      c_w(k,i) = 0.0_wp
[106]556                   ELSEIF ( c_w(k,i) > c_max )  THEN
557                      c_w(k,i) = c_max
558                   ENDIF
559                ELSE
560                   c_w(k,i) = c_max
[73]561                ENDIF
[106]562
[978]563                c_u_m_l(k) = c_u_m_l(k) + c_u(k,i)
564                c_v_m_l(k) = c_v_m_l(k) + c_v(k,i)
565                c_w_m_l(k) = c_w_m_l(k) + c_w(k,i)
[106]566
[978]567             ENDDO
568          ENDDO
[73]569
[978]570#if defined( __parallel )   
571          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
572          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
573                              MPI_SUM, comm1dx, ierr )   
574          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
575          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
576                              MPI_SUM, comm1dx, ierr ) 
577          IF ( collective_wait )  CALL MPI_BARRIER( comm1dx, ierr )
578          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
579                              MPI_SUM, comm1dx, ierr ) 
580#else
581          c_u_m = c_u_m_l
582          c_v_m = c_v_m_l
583          c_w_m = c_w_m_l
584#endif
585
586          c_u_m = c_u_m / (nx+1)
587          c_v_m = c_v_m / (nx+1)
588          c_w_m = c_w_m / (nx+1)
589
[73]590!
[978]591!--       Save old timelevels for the next timestep
592          IF ( intermediate_timestep_count == 1 )  THEN
593                u_m_n(:,:,:) = u(:,ny-1:ny,:)
594                v_m_n(:,:,:) = v(:,ny-1:ny,:)
595                w_m_n(:,:,:) = w(:,ny-1:ny,:)
596          ENDIF
[73]597
[978]598!
599!--       Calculate the new velocities
[996]600          DO  k = nzb+1, nzt+1
601             DO  i = nxlg, nxrg
[978]602                u_p(k,ny+1,i) = u(k,ny+1,i) - dt_3d * tsc(2) * c_u_m(k) *      &
603                                       ( u(k,ny+1,i) - u(k,ny,i) ) * ddy
[73]604
[978]605                v_p(k,ny+1,i) = v(k,ny+1,i)  - dt_3d * tsc(2) * c_v_m(k) *     &
606                                       ( v(k,ny+1,i) - v(k,ny,i) ) * ddy
[73]607
[978]608                w_p(k,ny+1,i) = w(k,ny+1,i) - dt_3d * tsc(2) * c_w_m(k) *      &
609                                       ( w(k,ny+1,i) - w(k,ny,i) ) * ddy
610             ENDDO
[1]611          ENDDO
612
613!
[978]614!--       Bottom boundary at the outflow
615          IF ( ibc_uv_b == 0 )  THEN
[1353]616             u_p(nzb,ny+1,:) = 0.0_wp
617             v_p(nzb,ny+1,:) = 0.0_wp   
[4268]618          ELSE
[978]619             u_p(nzb,ny+1,:) =  u_p(nzb+1,ny+1,:)
620             v_p(nzb,ny+1,:) =  v_p(nzb+1,ny+1,:)
621          ENDIF
[1353]622          w_p(nzb,ny+1,:) = 0.0_wp
[73]623
624!
[978]625!--       Top boundary at the outflow
626          IF ( ibc_uv_t == 0 )  THEN
627             u_p(nzt+1,ny+1,:) = u_init(nzt+1)
628             v_p(nzt+1,ny+1,:) = v_init(nzt+1)
629          ELSE
630             u_p(nzt+1,ny+1,:) = u_p(nzt,nyn+1,:)
631             v_p(nzt+1,ny+1,:) = v_p(nzt,nyn+1,:)
632          ENDIF
[1353]633          w_p(nzt:nzt+1,ny+1,:) = 0.0_wp
[978]634
[1]635       ENDIF
636
[75]637    ENDIF
638
[3182]639    IF ( bc_radiation_l )  THEN
[75]640
[1159]641       IF ( use_cmax )  THEN
[1717]642          u_p(:,:,0)  = u(:,:,1)
643          v_p(:,:,-1) = v(:,:,0)
[4268]644          w_p(:,:,-1) = w(:,:,0)
[1159]645       ELSEIF ( .NOT. use_cmax )  THEN
[75]646
[978]647          c_max = dx / dt_3d
[75]648
[1353]649          c_u_m_l = 0.0_wp 
650          c_v_m_l = 0.0_wp
651          c_w_m_l = 0.0_wp
[978]652
[1353]653          c_u_m = 0.0_wp 
654          c_v_m = 0.0_wp
655          c_w_m = 0.0_wp
[978]656
[1]657!
[996]658!--       Calculate the phase speeds for u, v, and w, first local and then
659!--       average along the outflow boundary.
660          DO  k = nzb+1, nzt+1
661             DO  j = nys, nyn
[75]662
[106]663                denom = u_m_l(k,j,1) - u_m_l(k,j,2)
664
[1353]665                IF ( denom /= 0.0_wp )  THEN
[996]666                   c_u(k,j) = -c_max * ( u(k,j,1) - u_m_l(k,j,1) ) / ( denom * tsc(2) )
[1353]667                   IF ( c_u(k,j) < 0.0_wp )  THEN
668                      c_u(k,j) = 0.0_wp
[107]669                   ELSEIF ( c_u(k,j) > c_max )  THEN
670                      c_u(k,j) = c_max
[106]671                   ENDIF
672                ELSE
[107]673                   c_u(k,j) = c_max
[75]674                ENDIF
675
[106]676                denom = v_m_l(k,j,0) - v_m_l(k,j,1)
[75]677
[1353]678                IF ( denom /= 0.0_wp )  THEN
[996]679                   c_v(k,j) = -c_max * ( v(k,j,0) - v_m_l(k,j,0) ) / ( denom * tsc(2) )
[1353]680                   IF ( c_v(k,j) < 0.0_wp )  THEN
681                      c_v(k,j) = 0.0_wp
[106]682                   ELSEIF ( c_v(k,j) > c_max )  THEN
683                      c_v(k,j) = c_max
684                   ENDIF
685                ELSE
686                   c_v(k,j) = c_max
[75]687                ENDIF
688
[106]689                denom = w_m_l(k,j,0) - w_m_l(k,j,1)
[75]690
[1353]691                IF ( denom /= 0.0_wp )  THEN
[996]692                   c_w(k,j) = -c_max * ( w(k,j,0) - w_m_l(k,j,0) ) / ( denom * tsc(2) )
[1353]693                   IF ( c_w(k,j) < 0.0_wp )  THEN
694                      c_w(k,j) = 0.0_wp
[106]695                   ELSEIF ( c_w(k,j) > c_max )  THEN
696                      c_w(k,j) = c_max
697                   ENDIF
698                ELSE
699                   c_w(k,j) = c_max
[75]700                ENDIF
[106]701
[978]702                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
703                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
704                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
[106]705
[978]706             ENDDO
707          ENDDO
[75]708
[978]709#if defined( __parallel )   
710          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
711          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
712                              MPI_SUM, comm1dy, ierr )   
713          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
714          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
715                              MPI_SUM, comm1dy, ierr ) 
716          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
717          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
718                              MPI_SUM, comm1dy, ierr ) 
719#else
720          c_u_m = c_u_m_l
721          c_v_m = c_v_m_l
722          c_w_m = c_w_m_l
723#endif
724
725          c_u_m = c_u_m / (ny+1)
726          c_v_m = c_v_m / (ny+1)
727          c_w_m = c_w_m / (ny+1)
728
[73]729!
[978]730!--       Save old timelevels for the next timestep
731          IF ( intermediate_timestep_count == 1 )  THEN
732                u_m_l(:,:,:) = u(:,:,1:2)
733                v_m_l(:,:,:) = v(:,:,0:1)
734                w_m_l(:,:,:) = w(:,:,0:1)
735          ENDIF
736
737!
738!--       Calculate the new velocities
[996]739          DO  k = nzb+1, nzt+1
[1113]740             DO  j = nysg, nyng
[978]741                u_p(k,j,0) = u(k,j,0) - dt_3d * tsc(2) * c_u_m(k) *            &
[106]742                                       ( u(k,j,0) - u(k,j,1) ) * ddx
[75]743
[978]744                v_p(k,j,-1) = v(k,j,-1) - dt_3d * tsc(2) * c_v_m(k) *          &
[75]745                                       ( v(k,j,-1) - v(k,j,0) ) * ddx
746
[978]747                w_p(k,j,-1) = w(k,j,-1) - dt_3d * tsc(2) * c_w_m(k) *          &
[75]748                                       ( w(k,j,-1) - w(k,j,0) ) * ddx
[978]749             ENDDO
[75]750          ENDDO
751
752!
[978]753!--       Bottom boundary at the outflow
754          IF ( ibc_uv_b == 0 )  THEN
[1353]755             u_p(nzb,:,0)  = 0.0_wp 
756             v_p(nzb,:,-1) = 0.0_wp
[4268]757          ELSE
[978]758             u_p(nzb,:,0)  =  u_p(nzb+1,:,0)
759             v_p(nzb,:,-1) =  v_p(nzb+1,:,-1)
760          ENDIF
[1353]761          w_p(nzb,:,-1) = 0.0_wp
[1]762
[75]763!
[978]764!--       Top boundary at the outflow
765          IF ( ibc_uv_t == 0 )  THEN
[1764]766             u_p(nzt+1,:,0)  = u_init(nzt+1)
[978]767             v_p(nzt+1,:,-1) = v_init(nzt+1)
768          ELSE
[1764]769             u_p(nzt+1,:,0)  = u_p(nzt,:,0)
[978]770             v_p(nzt+1,:,-1) = v_p(nzt,:,-1)
771          ENDIF
[1353]772          w_p(nzt:nzt+1,:,-1) = 0.0_wp
[978]773
[75]774       ENDIF
[73]775
[75]776    ENDIF
[73]777
[3182]778    IF ( bc_radiation_r )  THEN
[73]779
[1159]780       IF ( use_cmax )  THEN
781          u_p(:,:,nx+1) = u(:,:,nx)
782          v_p(:,:,nx+1) = v(:,:,nx)
[4268]783          w_p(:,:,nx+1) = w(:,:,nx)
[1159]784       ELSEIF ( .NOT. use_cmax )  THEN
[75]785
[978]786          c_max = dx / dt_3d
[75]787
[1353]788          c_u_m_l = 0.0_wp 
789          c_v_m_l = 0.0_wp
790          c_w_m_l = 0.0_wp
[978]791
[1353]792          c_u_m = 0.0_wp 
793          c_v_m = 0.0_wp
794          c_w_m = 0.0_wp
[978]795
[1]796!
[996]797!--       Calculate the phase speeds for u, v, and w, first local and then
798!--       average along the outflow boundary.
799          DO  k = nzb+1, nzt+1
800             DO  j = nys, nyn
[73]801
[106]802                denom = u_m_r(k,j,nx) - u_m_r(k,j,nx-1)
803
[1353]804                IF ( denom /= 0.0_wp )  THEN
[996]805                   c_u(k,j) = -c_max * ( u(k,j,nx) - u_m_r(k,j,nx) ) / ( denom * tsc(2) )
[1353]806                   IF ( c_u(k,j) < 0.0_wp )  THEN
807                      c_u(k,j) = 0.0_wp
[106]808                   ELSEIF ( c_u(k,j) > c_max )  THEN
809                      c_u(k,j) = c_max
810                   ENDIF
811                ELSE
812                   c_u(k,j) = c_max
[73]813                ENDIF
814
[106]815                denom = v_m_r(k,j,nx) - v_m_r(k,j,nx-1)
[73]816
[1353]817                IF ( denom /= 0.0_wp )  THEN
[996]818                   c_v(k,j) = -c_max * ( v(k,j,nx) - v_m_r(k,j,nx) ) / ( denom * tsc(2) )
[1353]819                   IF ( c_v(k,j) < 0.0_wp )  THEN
820                      c_v(k,j) = 0.0_wp
[106]821                   ELSEIF ( c_v(k,j) > c_max )  THEN
822                      c_v(k,j) = c_max
823                   ENDIF
824                ELSE
825                   c_v(k,j) = c_max
[73]826                ENDIF
827
[106]828                denom = w_m_r(k,j,nx) - w_m_r(k,j,nx-1)
[73]829
[1353]830                IF ( denom /= 0.0_wp )  THEN
[996]831                   c_w(k,j) = -c_max * ( w(k,j,nx) - w_m_r(k,j,nx) ) / ( denom * tsc(2) )
[1353]832                   IF ( c_w(k,j) < 0.0_wp )  THEN
833                      c_w(k,j) = 0.0_wp
[106]834                   ELSEIF ( c_w(k,j) > c_max )  THEN
835                      c_w(k,j) = c_max
836                   ENDIF
837                ELSE
838                   c_w(k,j) = c_max
[73]839                ENDIF
[106]840
[978]841                c_u_m_l(k) = c_u_m_l(k) + c_u(k,j)
842                c_v_m_l(k) = c_v_m_l(k) + c_v(k,j)
843                c_w_m_l(k) = c_w_m_l(k) + c_w(k,j)
[106]844
[978]845             ENDDO
846          ENDDO
[73]847
[4272]848#if defined( __parallel )
[978]849          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
850          CALL MPI_ALLREDUCE( c_u_m_l(nzb+1), c_u_m(nzb+1), nzt-nzb, MPI_REAL, &
[4272]851                              MPI_SUM, comm1dy, ierr )
[978]852          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
853          CALL MPI_ALLREDUCE( c_v_m_l(nzb+1), c_v_m(nzb+1), nzt-nzb, MPI_REAL, &
[4272]854                              MPI_SUM, comm1dy, ierr )
[978]855          IF ( collective_wait )  CALL MPI_BARRIER( comm1dy, ierr )
856          CALL MPI_ALLREDUCE( c_w_m_l(nzb+1), c_w_m(nzb+1), nzt-nzb, MPI_REAL, &
[4272]857                              MPI_SUM, comm1dy, ierr )
[978]858#else
859          c_u_m = c_u_m_l
860          c_v_m = c_v_m_l
861          c_w_m = c_w_m_l
862#endif
863
864          c_u_m = c_u_m / (ny+1)
865          c_v_m = c_v_m / (ny+1)
866          c_w_m = c_w_m / (ny+1)
867
[73]868!
[978]869!--       Save old timelevels for the next timestep
870          IF ( intermediate_timestep_count == 1 )  THEN
871                u_m_r(:,:,:) = u(:,:,nx-1:nx)
872                v_m_r(:,:,:) = v(:,:,nx-1:nx)
873                w_m_r(:,:,:) = w(:,:,nx-1:nx)
874          ENDIF
[73]875
[978]876!
877!--       Calculate the new velocities
[996]878          DO  k = nzb+1, nzt+1
[1113]879             DO  j = nysg, nyng
[978]880                u_p(k,j,nx+1) = u(k,j,nx+1) - dt_3d * tsc(2) * c_u_m(k) *      &
881                                       ( u(k,j,nx+1) - u(k,j,nx) ) * ddx
[73]882
[978]883                v_p(k,j,nx+1) = v(k,j,nx+1) - dt_3d * tsc(2) * c_v_m(k) *      &
884                                       ( v(k,j,nx+1) - v(k,j,nx) ) * ddx
[73]885
[978]886                w_p(k,j,nx+1) = w(k,j,nx+1) - dt_3d * tsc(2) * c_w_m(k) *      &
887                                       ( w(k,j,nx+1) - w(k,j,nx) ) * ddx
888             ENDDO
[73]889          ENDDO
890
891!
[978]892!--       Bottom boundary at the outflow
893          IF ( ibc_uv_b == 0 )  THEN
[1353]894             u_p(nzb,:,nx+1) = 0.0_wp
895             v_p(nzb,:,nx+1) = 0.0_wp 
[4268]896          ELSE
[978]897             u_p(nzb,:,nx+1) =  u_p(nzb+1,:,nx+1)
898             v_p(nzb,:,nx+1) =  v_p(nzb+1,:,nx+1)
899          ENDIF
[1353]900          w_p(nzb,:,nx+1) = 0.0_wp
[73]901
902!
[978]903!--       Top boundary at the outflow
904          IF ( ibc_uv_t == 0 )  THEN
905             u_p(nzt+1,:,nx+1) = u_init(nzt+1)
906             v_p(nzt+1,:,nx+1) = v_init(nzt+1)
907          ELSE
908             u_p(nzt+1,:,nx+1) = u_p(nzt,:,nx+1)
909             v_p(nzt+1,:,nx+1) = v_p(nzt,:,nx+1)
910          ENDIF
[1742]911          w_p(nzt:nzt+1,:,nx+1) = 0.0_wp
[978]912
[1]913       ENDIF
914
915    ENDIF
[3864]916
[1]917 END SUBROUTINE boundary_conds
Note: See TracBrowser for help on using the repository browser.