!-------------------------------------- LICENCE BEGIN ------------------------------------
!Environment Canada - Atmospheric Science and Technology License/Disclaimer, 
!                     version 3; Last Modified: May 7, 2008.
!This is free but copyrighted software; you can use/redistribute/modify it under the terms 
!of the Environment Canada - Atmospheric Science and Technology License/Disclaimer 
!version 3 or (at your option) any later version that should be found at: 
!http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html 
!
!This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; 
!without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. 
!See the above mentioned License/Disclaimer for more details.
!You should have received a copy of the License/Disclaimer along with this software; 
!if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec), 
!CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca
!-------------------------------------- LICENCE END --------------------------------------
***s/r vtopo_bnd_update - Update boundary nest values during orography change
*
#include "model_macros_f.h"
*

      subroutine vtopo_bnd_update(nest_q,nest_fi,nest_t,nest_pip, 1,1
     $  nest_s,DIST_DIM,nk,F_ndavx,F_ndavy)
*
      implicit none
      integer :: DIST_DIM,nk,F_ndavx,F_ndavy
      real, dimension(DIST_SHAPE) :: nest_s
      real, dimension(DIST_SHAPE,nk) :: nest_q,nest_fi,nest_t,nest_pip
*
*author
*     Ron McTaggart-Cowan - RPN - March 2007
*
*revision
* v3_30 - McTaggart-Cowan R.    - initial MPI version
*
*object
*     Determine the boundary status of the current MPI tile, and
*     call vtopo_predat to re-interpolate nesting fields in the
*     vertical in the blending zone.
*
*arguments
*  Name         I/O             Description
*--------------------------------------------------------------------------
* nk             I      Global z-dimension (G_nk)
* F_ndavx	 I	Thichness of sponge layer along x
* F_ndavy	 I	Thickness of sponge layer along y
* nest_q        I/O     Boundary conditions for log pressure (q=ln(p))
* nest_fi       I/O     Boundary conditions for phi (geopotential)
* nest_t        I/O     Boundary conditions for temperature
* nest_pip      I/O     Boundary conditions for perturbation pi
* nest_s        I/O     Boundary conditions for surface pi (ln (pi_s/z_s))
*--------------------------------------------------------------------------
*
*implicits
#include "glb_ld.cdk"


#include "out.cdk"


**
      type region
 	integer :: xmin,xmax,ymin,ymax
      end type
      integer, parameter :: un=1
      integer :: nit,njt,il,ih,jl,jh,is,js,nRegions,r
      type(region) :: n,e,s,w,ne,se,sw,nw
      type(region), dimension(8) :: zone
*
*     ---------------------------------------------------------------------
*
* If no blending zones are requested, return immediately
      if (F_ndavx == 0 .and. F_ndavy == 0) return
*
* No horizontal staggering in this implementation
      is = 0; js = 0
*
* Set blending zone extents
      nit = l_ni - is - pil_e
      njt = l_nj - js - pil_n
      il = 1 + F_ndavx + pil_w - 1
      ih = nit - F_ndavx + 1
      jl = 1 + F_ndavy + pil_s - 1
      jh = njt - F_ndavy + 1
*
* Define regions
      n%xmin = il; 		n%xmax = ih
      n%ymin = jh+un; 		n%ymax = njt
      e%xmin = ih+un; 		e%xmax = nit
      e%ymin = jl; 		e%ymax = jh
      s%xmin = il; 		s%xmax = ih
      s%ymin = 1+pil_s; 	s%ymax = jl-un
      w%xmin = 1+pil_w; 	w%xmax = il-un
      w%ymin = jl; 		w%ymax = jh
      ne%xmin = ih+un; 		ne%xmax = nit
      ne%ymin = jh+un; 		ne%ymax = njt
      se%xmin = ih+un; 		se%xmax = nit
      se%ymin = 1+pil_s; 	se%ymax = jl-un
      sw%xmin = 1+pil_w; 	sw%xmax = il-un
      sw%ymin = 1+pil_s; 	sw%ymax = jl-un
      nw%xmin = 1+pil_w; 	nw%xmax = il-un
      nw%ymin = jh+un; 		nw%ymax = njt
*
* Initialize the number of regions
      nRegions = -1
*
* Determine extents for all possible blending zones on the tile
      select case (north+south+east+west)
*
* All blending on this tile (4 zones)
      case(4)
*
*	Add blending zones to processing list
	nRegions = 8
	zone = (/n,e,s,w,ne,se,sw,nw/)
*
* Three blending zones on this tile
      case(3)
*
*	Add blending zones to processing list
	nRegions = 5; r = 1
	if (l_north .and. l_south) then
	  if (l_west) then
	    zone(r) = w; r=r+1
	    zone(r) = sw; r=r+1
	    zone(r) = nw; r=r+1
	    n%xmax = l_ni; s%xmax = l_ni
	  elseif (l_east) then
	    zone(r) = e; r=r+1
	    zone(r) = se; r=r+1
	    zone(r) = ne; r=r+1
	    n%xmin = 1; s%xmin = 1
	  endif
	  zone(r) = n; r=r+1
	  zone(r) = s; r=r+1
	else
	  if (l_north) then
	    zone(r) = n; r=r+1
	    zone(r) = ne; r=r+1
	    zone(r) = nw; r=r+1
	    e%ymin = 1; w%ymin = 1
	  elseif (l_south) then
	    zone(r) = s; r=r+1
	    zone(r) = se; r=r+1
	    zone(r) = sw; r=r+1
	    e%ymax = l_nj; w%ymax = l_nj
	  endif
	  zone(r) = e; r=r+1
	  zone(r) = w; r=r+1
	endif
*
* Two blending zones on this tile
      case(2)
*
*	Add blending zones to processing list
	r = 1
*       Check if we've got a corner piece
	if (((north+west) == 2).or.((north+east) == 2).or.
     $    ((south+west) == 2).or.((south+east) == 2)) then
	  nRegions = 3
	  if ((north+west) == 2) then
	    n%xmax = l_ni; w%ymin = 1
	    zone(r) = n; r=r+1
	    zone(r) = w; r=r+1
	    zone(r) = nw; r=r+1
	  elseif ((north+east) == 2) then
	    n%xmin = 1; e%ymin = 1
	    zone(r) = n; r=r+1
	    zone(r) = e; r=r+1
	    zone(r) = ne; r=r+1
	  elseif ((south+west) == 2) then
	    s%xmax = l_ni; w%ymax = l_nj
	    zone(r) = s; r=r+1
	    zone(r) = w; r=r+1
	    zone(r) = sw; r=r+1
	  elseif ((south+east) == 2) then
	    s%xmin = 1; e%ymax = l_nj
	    zone(r) = s; r=r+1
	    zone(r) = e; r=r+1
	    zone(r) = se; r=r+1
	  endif
*	No Corner piece - must be two opposite sides
	else
	  nRegions = 2
	  if (l_north) then
	    n%xmin = 1; n%xmax = l_ni
	    zone(r) = n; r=r+1
	  endif
	  if (l_east) then
	    e%ymin = 1; e%ymax = l_nj
	    zone(r) = e; r=r+1
	  endif
	  if (l_south) then
	    s%xmin = 1; s%xmax = l_ni
	    zone(r) = s; r=r+1
	  endif
	  if (l_west) then
	    w%ymin = 1; w%ymax = l_nj
	    zone(r) = w; r=r+1
	  endif
 	endif    
*
* One blending zone on this tile
      case(1)
*
*	Add blending zones to processing list
	nRegions = 1; r = 1
	if (l_north) then
	  n%xmin = 1; n%xmax = l_ni
	  zone(r) = n
	elseif (l_east) then
	  e%ymin = 1; e%ymax = l_nj
	  zone(r) = e
	elseif (l_south) then
	  s%xmin = 1; s%xmax = l_ni
	  zone(r) = s
	elseif (l_west) then
	  w%ymin = 1; w%ymax = l_nj
	  zone(r) = w
	endif
*
* No Blending zones on this tile
      end select
*
* Process regions
      do r=1,nRegions
	call vtopo_predat(nest_q,nest_fi,nest_t,nest_pip,nest_s,LDIST_DIM,
     $    zone(r)%xmin,zone(r)%xmax,zone(r)%ymin,zone(r)%ymax,nk)
      enddo
*
* End of subprogram
      return
      end subroutine vtopo_bnd_update