SUBROUTINE boundary_conds( range ) !------------------------------------------------------------------------------! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: boundary_conds.f90 768 2011-10-14 06:57:15Z franke $ ! ! 767 2011-10-14 06:39:12Z raasch ! ug,vg replaced by u_init,v_init as the Dirichlet top boundary condition ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng ! Removed mirror boundary conditions for u and v at the bottom in case of ! ibc_uv_b == 0. Instead, dirichelt boundary conditions (u=v=0) are set ! in init_3d_model ! ! 107 2007-08-17 13:54:45Z raasch ! Boundary conditions for temperature adjusted for coupled runs, ! bugfixes for the radiation boundary conditions at the outflow: radiation ! conditions are used for every substep, phase speeds are calculated for the ! first Runge-Kutta substep only and then reused, several index values changed ! ! 95 2007-06-02 16:48:38Z raasch ! Boundary conditions for salinity added ! ! 75 2007-03-22 09:54:05Z raasch ! The "main" part sets conditions for time level t+dt instead of level t, ! outflow boundary conditions changed from Neumann to radiation condition, ! uxrp, vynp eliminated, moisture renamed humidity ! ! 19 2007-02-23 04:53:48Z raasch ! Boundary conditions for e(nzt), pt(nzt), and q(nzt) removed because these ! gridpoints are now calculated by the prognostic equation, ! Dirichlet and zero gradient condition for pt established at top boundary ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.15 2006/02/23 09:54:55 raasch ! Surface boundary conditions in case of topography: nzb replaced by ! 2d-k-index-arrays (nzb_w_inner, etc.). Conditions for u and v remain ! unchanged (still using nzb) because a non-flat topography must use a ! Prandtl-layer, which don't requires explicit setting of the surface values. ! ! Revision 1.1 1997/09/12 06:21:34 raasch ! Initial revision ! ! ! Description: ! ------------ ! Boundary conditions for the prognostic quantities (range='main'). ! In case of non-cyclic lateral boundaries the conditions for velocities at ! the outflow are set after the pressure solver has been called (range= ! 'outflow_uvw'). ! One additional bottom boundary condition is applied for the TKE (=(u*)**2) ! in prandtl_fluxes. The cyclic lateral boundary conditions are implicitly ! handled in routine exchange_horiz. Pressure boundary conditions are ! explicitly set in routines pres, poisfft, poismg and sor. !------------------------------------------------------------------------------! USE arrays_3d USE control_parameters USE grid_variables USE indices USE pegrid IMPLICIT NONE CHARACTER (LEN=*) :: range INTEGER :: i, j, k REAL :: c_max, denom IF ( range == 'main') THEN ! !-- Bottom boundary IF ( ibc_uv_b == 1 ) THEN u_p(nzb,:,:) = u_p(nzb+1,:,:) v_p(nzb,:,:) = v_p(nzb+1,:,:) ENDIF DO i = nxlg, nxrg DO j = nysg, nyng w_p(nzb_w_inner(j,i),j,i) = 0.0 ENDDO ENDDO ! !-- Top boundary IF ( ibc_uv_t == 0 ) THEN u_p(nzt+1,:,:) = u_init(nzt+1) v_p(nzt+1,:,:) = v_init(nzt+1) ELSE u_p(nzt+1,:,:) = u_p(nzt,:,:) v_p(nzt+1,:,:) = v_p(nzt,:,:) ENDIF w_p(nzt:nzt+1,:,:) = 0.0 ! nzt is not a prognostic level (but cf. pres) ! !-- Temperature at bottom boundary. !-- In case of coupled runs (ibc_pt_b = 2) the temperature is given by !-- the sea surface temperature of the coupled ocean model. IF ( ibc_pt_b == 0 ) THEN DO i = nxlg, nxrg DO j = nysg, nyng pt_p(nzb_s_inner(j,i),j,i) = pt(nzb_s_inner(j,i),j,i) ENDDO ENDDO ELSEIF ( ibc_pt_b == 1 ) THEN DO i = nxlg, nxrg DO j = nysg, nyng pt_p(nzb_s_inner(j,i),j,i) = pt_p(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO ENDIF ! !-- Temperature at top boundary IF ( ibc_pt_t == 0 ) THEN pt_p(nzt+1,:,:) = pt(nzt+1,:,:) ELSEIF ( ibc_pt_t == 1 ) THEN pt_p(nzt+1,:,:) = pt_p(nzt,:,:) ELSEIF ( ibc_pt_t == 2 ) THEN pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1) ENDIF ! !-- Boundary conditions for TKE !-- Generally Neumann conditions with de/dz=0 are assumed IF ( .NOT. constant_diffusion ) THEN DO i = nxlg, nxrg DO j = nysg, nyng e_p(nzb_s_inner(j,i),j,i) = e_p(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO e_p(nzt+1,:,:) = e_p(nzt,:,:) ENDIF ! !-- Boundary conditions for salinity IF ( ocean ) THEN ! !-- Bottom boundary: Neumann condition because salinity flux is always !-- given DO i = nxlg, nxrg DO j = nysg, nyng sa_p(nzb_s_inner(j,i),j,i) = sa_p(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO ! !-- Top boundary: Dirichlet or Neumann IF ( ibc_sa_t == 0 ) THEN sa_p(nzt+1,:,:) = sa(nzt+1,:,:) ELSEIF ( ibc_sa_t == 1 ) THEN sa_p(nzt+1,:,:) = sa_p(nzt,:,:) ENDIF ENDIF ! !-- Boundary conditions for total water content or scalar, !-- bottom and top boundary (see also temperature) IF ( humidity .OR. passive_scalar ) THEN ! !-- Surface conditions for constant_humidity_flux IF ( ibc_q_b == 0 ) THEN DO i = nxlg, nxrg DO j = nysg, nyng q_p(nzb_s_inner(j,i),j,i) = q(nzb_s_inner(j,i),j,i) ENDDO ENDDO ELSE DO i = nxlg, nxrg DO j = nysg, nyng q_p(nzb_s_inner(j,i),j,i) = q_p(nzb_s_inner(j,i)+1,j,i) ENDDO ENDDO ENDIF ! !-- Top boundary q_p(nzt+1,:,:) = q_p(nzt,:,:) + bc_q_t_val * dzu(nzt+1) ENDIF ! !-- Lateral boundary conditions at the inflow. Quasi Neumann conditions !-- are needed for the wall normal velocity in order to ensure zero !-- divergence. Dirichlet conditions are used for all other quantities. IF ( inflow_s ) THEN v_p(:,nys,:) = v_p(:,nys-1,:) ELSEIF ( inflow_n ) THEN v_p(:,nyn,:) = v_p(:,nyn+1,:) ELSEIF ( inflow_l ) THEN u_p(:,:,nxl) = u_p(:,:,nxl-1) ELSEIF ( inflow_r ) THEN u_p(:,:,nxr) = u_p(:,:,nxr+1) ENDIF ! !-- Lateral boundary conditions for scalar quantities at the outflow IF ( outflow_s ) THEN pt_p(:,nys-1,:) = pt_p(:,nys,:) IF ( .NOT. constant_diffusion ) e_p(:,nys-1,:) = e_p(:,nys,:) IF ( humidity .OR. passive_scalar ) q_p(:,nys-1,:) = q_p(:,nys,:) ELSEIF ( outflow_n ) THEN pt_p(:,nyn+1,:) = pt_p(:,nyn,:) IF ( .NOT. constant_diffusion ) e_p(:,nyn+1,:) = e_p(:,nyn,:) IF ( humidity .OR. passive_scalar ) q_p(:,nyn+1,:) = q_p(:,nyn,:) ELSEIF ( outflow_l ) THEN pt_p(:,:,nxl-1) = pt_p(:,:,nxl) IF ( .NOT. constant_diffusion ) e_p(:,:,nxl-1) = e_p(:,:,nxl) IF ( humidity .OR. passive_scalar ) q_p(:,:,nxl-1) = q_p(:,:,nxl) ELSEIF ( outflow_r ) THEN pt_p(:,:,nxr+1) = pt_p(:,:,nxr) IF ( .NOT. constant_diffusion ) e_p(:,:,nxr+1) = e_p(:,:,nxr) IF ( humidity .OR. passive_scalar ) q_p(:,:,nxr+1) = q_p(:,:,nxr) ENDIF ENDIF ! !-- Radiation boundary condition for the velocities at the respective outflow IF ( outflow_s ) THEN c_max = dy / dt_3d DO i = nxlg, nxrg DO k = nzb+1, nzt+1 ! !-- Calculate the phase speeds for u,v, and w. In case of using a !-- Runge-Kutta scheme, do this for the first substep only and then !-- reuse this values for the further substeps. IF ( intermediate_timestep_count == 1 ) THEN denom = u_m_s(k,0,i) - u_m_s(k,1,i) IF ( denom /= 0.0 ) THEN c_u(k,i) = -c_max * ( u(k,0,i) - u_m_s(k,0,i) ) / denom IF ( c_u(k,i) < 0.0 ) THEN c_u(k,i) = 0.0 ELSEIF ( c_u(k,i) > c_max ) THEN c_u(k,i) = c_max ENDIF ELSE c_u(k,i) = c_max ENDIF denom = v_m_s(k,1,i) - v_m_s(k,2,i) IF ( denom /= 0.0 ) THEN c_v(k,i) = -c_max * ( v(k,1,i) - v_m_s(k,1,i) ) / denom IF ( c_v(k,i) < 0.0 ) THEN c_v(k,i) = 0.0 ELSEIF ( c_v(k,i) > c_max ) THEN c_v(k,i) = c_max ENDIF ELSE c_v(k,i) = c_max ENDIF denom = w_m_s(k,0,i) - w_m_s(k,1,i) IF ( denom /= 0.0 ) THEN c_w(k,i) = -c_max * ( w(k,0,i) - w_m_s(k,0,i) ) / denom IF ( c_w(k,i) < 0.0 ) THEN c_w(k,i) = 0.0 ELSEIF ( c_w(k,i) > c_max ) THEN c_w(k,i) = c_max ENDIF ELSE c_w(k,i) = c_max ENDIF ! !-- Save old timelevels for the next timestep u_m_s(k,:,i) = u(k,0:1,i) v_m_s(k,:,i) = v(k,1:2,i) w_m_s(k,:,i) = w(k,0:1,i) ENDIF ! !-- Calculate the new velocities u_p(k,-1,i) = u(k,-1,i) - dt_3d * tsc(2) * c_u(k,i) * & ( u(k,-1,i) - u(k,0,i) ) * ddy v_p(k,0,i) = v(k,0,i) - dt_3d * tsc(2) * c_v(k,i) * & ( v(k,0,i) - v(k,1,i) ) * ddy w_p(k,-1,i) = w(k,-1,i) - dt_3d * tsc(2) * c_w(k,i) * & ( w(k,-1,i) - w(k,0,i) ) * ddy ENDDO ENDDO ! !-- Bottom boundary at the outflow IF ( ibc_uv_b == 0 ) THEN u_p(nzb,-1,:) = 0.0 v_p(nzb,0,:) = 0.0 ELSE u_p(nzb,-1,:) = u_p(nzb+1,-1,:) v_p(nzb,0,:) = v_p(nzb+1,0,:) ENDIF w_p(nzb,-1,:) = 0.0 ! !-- Top boundary at the outflow IF ( ibc_uv_t == 0 ) THEN u_p(nzt+1,-1,:) = u_init(nzt+1) v_p(nzt+1,0,:) = v_init(nzt+1) ELSE u_p(nzt+1,-1,:) = u(nzt,-1,:) v_p(nzt+1,0,:) = v(nzt,0,:) ENDIF w_p(nzt:nzt+1,-1,:) = 0.0 ENDIF IF ( outflow_n ) THEN c_max = dy / dt_3d DO i = nxlg, nxrg DO k = nzb+1, nzt+1 ! !-- Calculate the phase speeds for u,v, and w. In case of using a !-- Runge-Kutta scheme, do this for the first substep only and then !-- reuse this values for the further substeps. IF ( intermediate_timestep_count == 1 ) THEN denom = u_m_n(k,ny,i) - u_m_n(k,ny-1,i) IF ( denom /= 0.0 ) THEN c_u(k,i) = -c_max * ( u(k,ny,i) - u_m_n(k,ny,i) ) / denom IF ( c_u(k,i) < 0.0 ) THEN c_u(k,i) = 0.0 ELSEIF ( c_u(k,i) > c_max ) THEN c_u(k,i) = c_max ENDIF ELSE c_u(k,i) = c_max ENDIF denom = v_m_n(k,ny,i) - v_m_n(k,ny-1,i) IF ( denom /= 0.0 ) THEN c_v(k,i) = -c_max * ( v(k,ny,i) - v_m_n(k,ny,i) ) / denom IF ( c_v(k,i) < 0.0 ) THEN c_v(k,i) = 0.0 ELSEIF ( c_v(k,i) > c_max ) THEN c_v(k,i) = c_max ENDIF ELSE c_v(k,i) = c_max ENDIF denom = w_m_n(k,ny,i) - w_m_n(k,ny-1,i) IF ( denom /= 0.0 ) THEN c_w(k,i) = -c_max * ( w(k,ny,i) - w_m_n(k,ny,i) ) / denom IF ( c_w(k,i) < 0.0 ) THEN c_w(k,i) = 0.0 ELSEIF ( c_w(k,i) > c_max ) THEN c_w(k,i) = c_max ENDIF ELSE c_w(k,i) = c_max ENDIF ! !-- Swap timelevels for the next timestep u_m_n(k,:,i) = u(k,ny-1:ny,i) v_m_n(k,:,i) = v(k,ny-1:ny,i) w_m_n(k,:,i) = w(k,ny-1:ny,i) ENDIF ! !-- Calculate the new velocities u_p(k,ny+1,i) = u(k,ny+1,i) - dt_3d * tsc(2) * c_u(k,i) * & ( u(k,ny+1,i) - u(k,ny,i) ) * ddy v_p(k,ny+1,i) = v(k,ny+1,i) - dt_3d * tsc(2) * c_v(k,i) * & ( v(k,ny+1,i) - v(k,ny,i) ) * ddy w_p(k,ny+1,i) = w(k,ny+1,i) - dt_3d * tsc(2) * c_w(k,i) * & ( w(k,ny+1,i) - w(k,ny,i) ) * ddy ENDDO ENDDO ! !-- Bottom boundary at the outflow IF ( ibc_uv_b == 0 ) THEN u_p(nzb,ny+1,:) = 0.0 v_p(nzb,ny+1,:) = 0.0 ELSE u_p(nzb,ny+1,:) = u_p(nzb+1,ny+1,:) v_p(nzb,ny+1,:) = v_p(nzb+1,ny+1,:) ENDIF w_p(nzb,ny+1,:) = 0.0 ! !-- Top boundary at the outflow IF ( ibc_uv_t == 0 ) THEN u_p(nzt+1,ny+1,:) = u_init(nzt+1) v_p(nzt+1,ny+1,:) = v_init(nzt+1) ELSE u_p(nzt+1,ny+1,:) = u_p(nzt,nyn+1,:) v_p(nzt+1,ny+1,:) = v_p(nzt,nyn+1,:) ENDIF w_p(nzt:nzt+1,ny+1,:) = 0.0 ENDIF IF ( outflow_l ) THEN c_max = dx / dt_3d DO j = nysg, nyng DO k = nzb+1, nzt+1 ! !-- Calculate the phase speeds for u,v, and w. In case of using a !-- Runge-Kutta scheme, do this for the first substep only and then !-- reuse this values for the further substeps. IF ( intermediate_timestep_count == 1 ) THEN denom = u_m_l(k,j,1) - u_m_l(k,j,2) IF ( denom /= 0.0 ) THEN c_u(k,j) = -c_max * ( u(k,j,1) - u_m_l(k,j,1) ) / denom IF ( c_u(k,j) < 0.0 ) THEN c_u(k,j) = 0.0 ELSEIF ( c_u(k,j) > c_max ) THEN c_u(k,j) = c_max ENDIF ELSE c_u(k,j) = c_max ENDIF denom = v_m_l(k,j,0) - v_m_l(k,j,1) IF ( denom /= 0.0 ) THEN c_v(k,j) = -c_max * ( v(k,j,0) - v_m_l(k,j,0) ) / denom IF ( c_v(k,j) < 0.0 ) THEN c_v(k,j) = 0.0 ELSEIF ( c_v(k,j) > c_max ) THEN c_v(k,j) = c_max ENDIF ELSE c_v(k,j) = c_max ENDIF denom = w_m_l(k,j,0) - w_m_l(k,j,1) IF ( denom /= 0.0 ) THEN c_w(k,j) = -c_max * ( w(k,j,0) - w_m_l(k,j,0) ) / denom IF ( c_w(k,j) < 0.0 ) THEN c_w(k,j) = 0.0 ELSEIF ( c_w(k,j) > c_max ) THEN c_w(k,j) = c_max ENDIF ELSE c_w(k,j) = c_max ENDIF ! !-- Swap timelevels for the next timestep u_m_l(k,j,:) = u(k,j,1:2) v_m_l(k,j,:) = v(k,j,0:1) w_m_l(k,j,:) = w(k,j,0:1) ENDIF ! !-- Calculate the new velocities u_p(k,j,0) = u(k,j,0) - dt_3d * tsc(2) * c_u(k,j) * & ( u(k,j,0) - u(k,j,1) ) * ddx v_p(k,j,-1) = v(k,j,-1) - dt_3d * tsc(2) * c_v(k,j) * & ( v(k,j,-1) - v(k,j,0) ) * ddx w_p(k,j,-1) = w(k,j,-1) - dt_3d * tsc(2) * c_w(k,j) * & ( w(k,j,-1) - w(k,j,0) ) * ddx ENDDO ENDDO ! !-- Bottom boundary at the outflow IF ( ibc_uv_b == 0 ) THEN u_p(nzb,:,0) = 0.0 v_p(nzb,:,-1) = 0.0 ELSE u_p(nzb,:,0) = u_p(nzb+1,:,0) v_p(nzb,:,-1) = v_p(nzb+1,:,-1) ENDIF w_p(nzb,:,-1) = 0.0 ! !-- Top boundary at the outflow IF ( ibc_uv_t == 0 ) THEN u_p(nzt+1,:,-1) = u_init(nzt+1) v_p(nzt+1,:,-1) = v_init(nzt+1) ELSE u_p(nzt+1,:,-1) = u_p(nzt,:,-1) v_p(nzt+1,:,-1) = v_p(nzt,:,-1) ENDIF w_p(nzt:nzt+1,:,-1) = 0.0 ENDIF IF ( outflow_r ) THEN c_max = dx / dt_3d DO j = nysg, nyng DO k = nzb+1, nzt+1 ! !-- Calculate the phase speeds for u,v, and w. In case of using a !-- Runge-Kutta scheme, do this for the first substep only and then !-- reuse this values for the further substeps. IF ( intermediate_timestep_count == 1 ) THEN denom = u_m_r(k,j,nx) - u_m_r(k,j,nx-1) IF ( denom /= 0.0 ) THEN c_u(k,j) = -c_max * ( u(k,j,nx) - u_m_r(k,j,nx) ) / denom IF ( c_u(k,j) < 0.0 ) THEN c_u(k,j) = 0.0 ELSEIF ( c_u(k,j) > c_max ) THEN c_u(k,j) = c_max ENDIF ELSE c_u(k,j) = c_max ENDIF denom = v_m_r(k,j,nx) - v_m_r(k,j,nx-1) IF ( denom /= 0.0 ) THEN c_v(k,j) = -c_max * ( v(k,j,nx) - v_m_r(k,j,nx) ) / denom IF ( c_v(k,j) < 0.0 ) THEN c_v(k,j) = 0.0 ELSEIF ( c_v(k,j) > c_max ) THEN c_v(k,j) = c_max ENDIF ELSE c_v(k,j) = c_max ENDIF denom = w_m_r(k,j,nx) - w_m_r(k,j,nx-1) IF ( denom /= 0.0 ) THEN c_w(k,j) = -c_max * ( w(k,j,nx) - w_m_r(k,j,nx) ) / denom IF ( c_w(k,j) < 0.0 ) THEN c_w(k,j) = 0.0 ELSEIF ( c_w(k,j) > c_max ) THEN c_w(k,j) = c_max ENDIF ELSE c_w(k,j) = c_max ENDIF ! !-- Swap timelevels for the next timestep u_m_r(k,j,:) = u(k,j,nx-1:nx) v_m_r(k,j,:) = v(k,j,nx-1:nx) w_m_r(k,j,:) = w(k,j,nx-1:nx) ENDIF ! !-- Calculate the new velocities u_p(k,j,nx+1) = u(k,j,nx+1) - dt_3d * tsc(2) * c_u(k,j) * & ( u(k,j,nx+1) - u(k,j,nx) ) * ddx v_p(k,j,nx+1) = v(k,j,nx+1) - dt_3d * tsc(2) * c_v(k,j) * & ( v(k,j,nx+1) - v(k,j,nx) ) * ddx w_p(k,j,nx+1) = w(k,j,nx+1) - dt_3d * tsc(2) * c_w(k,j) * & ( w(k,j,nx+1) - w(k,j,nx) ) * ddx ENDDO ENDDO ! !-- Bottom boundary at the outflow IF ( ibc_uv_b == 0 ) THEN u_p(nzb,:,nx+1) = 0.0 v_p(nzb,:,nx+1) = 0.0 ELSE u_p(nzb,:,nx+1) = u_p(nzb+1,:,nx+1) v_p(nzb,:,nx+1) = v_p(nzb+1,:,nx+1) ENDIF w_p(nzb,:,nx+1) = 0.0 ! !-- Top boundary at the outflow IF ( ibc_uv_t == 0 ) THEN u_p(nzt+1,:,nx+1) = u_init(nzt+1) v_p(nzt+1,:,nx+1) = v_init(nzt+1) ELSE u_p(nzt+1,:,nx+1) = u_p(nzt,:,nx+1) v_p(nzt+1,:,nx+1) = v_p(nzt,:,nx+1) ENDIF w(nzt:nzt+1,:,nx+1) = 0.0 ENDIF END SUBROUTINE boundary_conds