cpld_gridgen  1.13.0
 All Data Structures Files Functions Variables Pages
vertices.F90
Go to the documentation of this file.
1 
12 
13 module vertices
14 
15  use gengrid_kinds, only : dbl_kind
16  use grdvars, only : ni,nj,nv
17 
18  implicit none
19 
20 contains
32 
33  subroutine fill_vertices(jbeg,jend,iVert,jVert,lat,lon,latvert,lonvert)
34 
35  integer, intent( in) :: jbeg,jend
36  integer, intent( in) :: ivert(nv), jvert(nv)
37  real(dbl_kind), dimension(ni,nj), intent( in) :: lat, lon
38 
39  real(dbl_kind), dimension(ni,nj,nv), intent(out) :: latvert,lonvert
40 
41  ! local variables
42  integer :: i,j,n,ii,jj
43 
44  do j = jbeg,jend
45  do i = 1,ni
46  do n = 1,nv
47  ii = i + ivert(n); jj = j + jvert(n)
48  if(ii .eq. 0)ii = ni
49  if(ii .eq. ni+1)ii = 1
50  latvert(i,j,n) = lat(ii,jj)
51  lonvert(i,j,n) = lon(ii,jj)
52  enddo
53  enddo
54  enddo
55  end subroutine fill_vertices
56 
67 
68  subroutine fill_bottom(iVert,jVert,lat,lon,latvert,lonvert,dlat)
69 
70  integer, intent( in) :: ivert(nv), jvert(nv)
71  real(dbl_kind), dimension(ni,nj), intent( in) :: lat, lon
72  real(dbl_kind), dimension(ni), intent( in) :: dlat
73 
74  real(dbl_kind), dimension(ni,nj,nv), intent(out) :: latvert,lonvert
75 
76  ! local variables
77  integer :: i,j,n,ii,jj
78 
79  ! fill in grid bottom (j=1)
80  ! vertices 1,2 are available
81  ! vertices 3,4 must be set manually
82  j = 1
83  do i = 1,ni
84  do n = 1,2
85  ii = i + ivert(n); jj = j + jvert(n)
86  if(ii .eq. 0)ii = ni
87  if(ii .eq. ni+1)ii = 1
88  latvert(i,j,n) = lat(ii,jj)
89  lonvert(i,j,n) = lon(ii,jj)
90  enddo
91  do n = 3,4
92  ii = i + ivert(n)
93  if(ii .eq. 0)ii = ni
94  if(ii .eq. ni+1)ii = 1
95  latvert(i,j, n) = dlat(ii)
96  enddo
97  lonvert(i,j, 3) = lonvert(i,j,2)
98  lonvert(i,j, 4) = lonvert(i,j,1)
99  enddo
100  end subroutine fill_bottom
101 
113 
114  subroutine fill_top(iVert,jVert,lat,lon,latvert,lonvert,xlat,xlon)
115 
116  integer, intent( in) :: ivert(nv), jvert(nv)
117  real(dbl_kind), dimension(ni,nj), intent( in) :: lat, lon
118  real(dbl_kind), dimension(ni), intent( in) :: xlat, xlon
119 
120  real(dbl_kind), dimension(ni,nj,nv), intent(out) :: latvert,lonvert
121 
122  ! local variables
123  integer :: i,j,n,ii,jj
124 
125  ! fill in grid top (j=nj)
126  ! vertices 3,4 are available
127  ! vertices 1,2 must be set manually using 'across seam' values
128  j = nj
129  do i = 1,ni
130  do n = 3,4
131  ii = i + ivert(n); jj = j + jvert(n)
132  if(ii .eq. 0)ii = ni
133  if(ii .eq. ni+1)ii = 1
134  latvert(i,j,n) = lat(ii,jj)
135  lonvert(i,j,n) = lon(ii,jj)
136  enddo
137  do n = 1,2
138  ii = i + ivert(n)
139  if(ii .eq. 0)ii = ni
140  if(ii .eq. ni+1)ii = 1
141  latvert(i,j,n) = xlat(ii)
142  lonvert(i,j,n) = xlon(ii)
143  enddo
144  enddo
145  !latCv_vert(i,j, 1) = latCv_vert(i,j,4)
146  !latCv_vert(i,j, 2) = latCv_vert(i,j,3)
147  !lonCv_vert(i,j, 1) = lonCv_vert(i,j,4)+240.d0
148  !lonCv_vert(i,j, 2) = lonCv_vert(i,j,3)+240.d0
149  end subroutine fill_top
150 end module vertices
subroutine fill_vertices(jbeg, jend, iVert, jVert, lat, lon, latvert, lonvert)
Fill the vertices for any stagger location between bounding j-rows.
Definition: vertices.F90:33
subroutine fill_bottom(iVert, jVert, lat, lon, latvert, lonvert, dlat)
Fill the vertices for a stagger location along the bottom j-row.
Definition: vertices.F90:68
subroutine fill_top(iVert, jVert, lat, lon, latvert, lonvert, xlat, xlon)
Fill the vertices for a stagger location along the top j-row.
Definition: vertices.F90:114