global_cycle 1.14.0
Loading...
Searching...
No Matches
utils.F90
Go to the documentation of this file.
1
4
7 module utils
8
9 public remap_coef
10
11 contains
12
38 SUBROUTINE remap_coef( is, ie, js, je,&
39 im, jm, lon, lat, id1, id2, jdc, s2c, agrid )
40
41 implicit none
42 integer, intent(in):: is, ie, js, je
43 integer, intent(in):: im, jm
44 real, intent(in):: lon(im), lat(jm)
45 real, intent(out):: s2c(is:ie,js:je,4)
46 integer, intent(out), dimension(is:ie,js:je):: id1, id2, jdc
47 real, intent(in):: agrid(is:ie,js:je,2)
48 ! local:
49 real :: rdlon(im)
50 real :: rdlat(jm)
51 real:: a1, b1
52 real, parameter :: pi = 3.1415926
53 integer i,j, i1, i2, jc, i0, j0
54 do i=1,im-1
55 rdlon(i) = 1. / (lon(i+1) - lon(i))
56 enddo
57 rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
58
59 do j=1,jm-1
60 rdlat(j) = 1. / (lat(j+1) - lat(j))
61 enddo
62
63 ! * Interpolate to cubed sphere cell center
64 do 5000 j=js,je
65
66 do i=is,ie
67
68 if ( agrid(i,j,1)>lon(im) ) then
69 i1 = im; i2 = 1
70 a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
71 elseif ( agrid(i,j,1)<lon(1) ) then
72 i1 = im; i2 = 1
73 a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
74 else
75 do i0=1,im-1
76 if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) ) then
77 i1 = i0; i2 = i0+1
78 a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
79 go to 111
80 endif
81 enddo
82 endif
83111 continue
84
85 if ( agrid(i,j,2)<lat(1) ) then
86 jc = 1
87 b1 = 0.
88 elseif ( agrid(i,j,2)>lat(jm) ) then
89 jc = jm-1
90 b1 = 1.
91 else
92 do j0=1,jm-1
93 if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) ) then
94 jc = j0
95 b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
96 go to 222
97 endif
98 enddo
99 endif
100222 continue
101
102 if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 ) then
103 write(*,*) 'gid=', i,j,a1, b1
104 endif
105
106 s2c(i,j,1) = (1.-a1) * (1.-b1)
107 s2c(i,j,2) = a1 * (1.-b1)
108 s2c(i,j,3) = a1 * b1
109 s2c(i,j,4) = (1.-a1) * b1
110 id1(i,j) = i1
111 id2(i,j) = i2
112 jdc(i,j) = jc
113 enddo !i-loop
1145000 continue ! j-loop
115
116 END SUBROUTINE remap_coef
117
118 end module utils
Module containing utility routines.
Definition utils.F90:7
subroutine, public remap_coef(is, ie, js, je, im, jm, lon, lat, id1, id2, jdc, s2c, agrid)
Generate the weights and index of the grids used in the bilinear interpolation.
Definition utils.F90:40