SCALE-RM
scale_intepolation.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14  !-----------------------------------------------------------------------------
15  !
16  !++ used modules
17  !
18  use scale_precision
19  use scale_stdio
20  use scale_prof
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: interp_setup
30  public :: interp_vertical_xi2z
31  public :: interp_vertical_z2xi
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  logical, public :: interp_available = .false.
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private procedure
42  !
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  integer, private, allocatable :: interp_xi2z_idx (:,:,:,:)
48  real(RP), private, allocatable :: interp_xi2z_coef(:,:,:,:)
49  integer, private, allocatable :: interp_z2xi_idx (:,:,:,:)
50  real(RP), private, allocatable :: interp_z2xi_coef(:,:,:,:)
51 
52  !-----------------------------------------------------------------------------
53 contains
54  !-----------------------------------------------------------------------------
56  subroutine interp_setup
57  use scale_grid, only: &
58  grid_cz, &
59  grid_fz
60  use scale_topography, only: &
62  use scale_grid_real, only: &
63  real_cz, &
64  real_fz
65  implicit none
66 
67  integer :: k, i, j, kk, kp
68  !---------------------------------------------------------------------------
69 
70  if( io_l ) write(io_fid_log,*)
71  if( io_l ) write(io_fid_log,*) '++++++ Module[INTERPOLATION] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
72  if( io_l ) write(io_fid_log,*) '*** No namelists.'
73 
75  if( io_l ) write(io_fid_log,*)
76  if( io_l ) write(io_fid_log,*) '*** Vertical interpolation is available? : ', interp_available
77 
78  allocate( interp_xi2z_idx(ka,ia,ja,2) )
79  allocate( interp_xi2z_coef(ka,ia,ja,3) )
80  allocate( interp_z2xi_idx(ka,ia,ja,2) )
81  allocate( interp_z2xi_coef(ka,ia,ja,3) )
82 
83  do j = 1, ja
84  do i = 1, ia
85  do k = ks, ke
86  if ( grid_cz(k) <= real_fz(ks-1,i,j) ) then
87 
88  interp_xi2z_idx(k,i,j,1) = ks ! dummmy
89  interp_xi2z_idx(k,i,j,2) = ks ! dummmy
90  interp_xi2z_coef(k,i,j,1) = 0.0_rp
91  interp_xi2z_coef(k,i,j,2) = 0.0_rp
92  interp_xi2z_coef(k,i,j,3) = 1.0_rp ! set UNDEF
93 
94  elseif( grid_cz(k) <= real_cz(ks,i,j) ) then
95 
96  interp_xi2z_idx(k,i,j,1) = ks ! dummmy
97  interp_xi2z_idx(k,i,j,2) = ks
98  interp_xi2z_coef(k,i,j,1) = 0.0_rp
99  interp_xi2z_coef(k,i,j,2) = 1.0_rp
100  interp_xi2z_coef(k,i,j,3) = 0.0_rp
101 
102  elseif( grid_cz(k) > real_cz(ke,i,j) ) then
103 
104  interp_xi2z_idx(k,i,j,1) = ke
105  interp_xi2z_idx(k,i,j,2) = ke ! dummmy
106  interp_xi2z_coef(k,i,j,1) = 1.0_rp
107  interp_xi2z_coef(k,i,j,2) = 0.0_rp
108  interp_xi2z_coef(k,i,j,3) = 0.0_rp
109 
110  elseif( grid_cz(k) > real_fz(ke,i,j) ) then
111 
112  interp_xi2z_idx(k,i,j,1) = ke ! dummmy
113  interp_xi2z_idx(k,i,j,2) = ke ! dummmy
114  interp_xi2z_coef(k,i,j,1) = 0.0_rp
115  interp_xi2z_coef(k,i,j,2) = 0.0_rp
116  interp_xi2z_coef(k,i,j,3) = 1.0_rp ! set UNDEF
117 
118  else
119 
120  do kk = ks+1, ke
121  kp = kk
122  if( grid_cz(k) <= real_cz(kk,i,j) ) exit
123  enddo
124 
125  interp_xi2z_idx(k,i,j,1) = kp - 1
126  interp_xi2z_idx(k,i,j,2) = kp
127  interp_xi2z_coef(k,i,j,1) = ( real_cz(kp,i,j) - grid_cz(k) ) &
128  / ( real_cz(kp,i,j) - real_cz(kp-1,i,j) )
129  interp_xi2z_coef(k,i,j,2) = ( grid_cz(k) - real_cz(kp-1,i,j) ) &
130  / ( real_cz(kp,i,j) - real_cz(kp-1,i,j) )
131  interp_xi2z_coef(k,i,j,3) = 0.0_rp
132 
133  endif
134  enddo
135  enddo
136  enddo
137 
138  do j = 1, ja
139  do i = 1, ia
140  do k = ks, ke
141  if ( real_cz(k,i,j) <= grid_fz(ks-1) ) then
142 
143  interp_z2xi_idx(k,i,j,1) = ks ! dummmy
144  interp_z2xi_idx(k,i,j,2) = ks ! dummmy
145  interp_z2xi_coef(k,i,j,1) = 0.0_rp
146  interp_z2xi_coef(k,i,j,2) = 0.0_rp
147  interp_z2xi_coef(k,i,j,3) = 1.0_rp ! set UNDEF
148 
149  elseif( real_cz(k,i,j) <= grid_cz(ks) ) then
150 
151  interp_z2xi_idx(k,i,j,1) = ks ! dummmy
152  interp_z2xi_idx(k,i,j,2) = ks
153  interp_z2xi_coef(k,i,j,1) = 0.0_rp
154  interp_z2xi_coef(k,i,j,2) = 1.0_rp
155  interp_z2xi_coef(k,i,j,3) = 0.0_rp
156 
157  elseif( real_cz(k,i,j) > grid_cz(ke) ) then
158 
159  interp_z2xi_idx(k,i,j,1) = ke
160  interp_z2xi_idx(k,i,j,2) = ke ! dummmy
161  interp_z2xi_coef(k,i,j,1) = 1.0_rp
162  interp_z2xi_coef(k,i,j,2) = 0.0_rp
163  interp_z2xi_coef(k,i,j,3) = 0.0_rp
164 
165  elseif( real_cz(k,i,j) > grid_fz(ke) ) then
166 
167  interp_z2xi_idx(k,i,j,1) = ke ! dummmy
168  interp_z2xi_idx(k,i,j,2) = ke ! dummmy
169  interp_z2xi_coef(k,i,j,1) = 0.0_rp
170  interp_z2xi_coef(k,i,j,2) = 0.0_rp
171  interp_z2xi_coef(k,i,j,3) = 1.0_rp ! set UNDEF
172 
173  else
174 
175  do kk = ks+1, ke
176  kp = kk
177  if( real_cz(k,i,j) <= grid_cz(kk) ) exit
178  enddo
179 
180  interp_z2xi_idx(k,i,j,1) = kp - 1
181  interp_z2xi_idx(k,i,j,2) = kp
182  interp_z2xi_coef(k,i,j,1) = ( grid_cz(kp) - real_cz(k,i,j) ) &
183  / ( grid_cz(kp) - grid_cz(kp-1) )
184  interp_z2xi_coef(k,i,j,2) = ( real_cz(k,i,j) - grid_cz(kp-1) ) &
185  / ( grid_cz(kp) - grid_cz(kp-1) )
186  interp_z2xi_coef(k,i,j,3) = 0.0_rp
187 
188  endif
189  enddo
190  enddo
191  enddo
192 
193  do j = 1, ja
194  do i = 1, ia
195  interp_xi2z_idx( 1:ks-1,i,j,1) = ks ! dummmy
196  interp_xi2z_idx( 1:ks-1,i,j,2) = ks ! dummmy
197  interp_xi2z_coef( 1:ks-1,i,j,1) = 0.0_rp
198  interp_xi2z_coef( 1:ks-1,i,j,2) = 0.0_rp
199  interp_xi2z_coef( 1:ks-1,i,j,3) = 1.0_rp ! set UNDEF
200 
201  interp_xi2z_idx(ke+1:ka,i,j,1) = ke ! dummmy
202  interp_xi2z_idx(ke+1:ka,i,j,2) = ke ! dummmy
203  interp_xi2z_coef(ke+1:ka,i,j,1) = 0.0_rp
204  interp_xi2z_coef(ke+1:ka,i,j,2) = 0.0_rp
205  interp_xi2z_coef(ke+1:ka,i,j,3) = 1.0_rp ! set UNDEF
206 
207  interp_z2xi_idx( 1:ks-1,i,j,1) = ks ! dummmy
208  interp_z2xi_idx( 1:ks-1,i,j,2) = ks ! dummmy
209  interp_z2xi_coef( 1:ks-1,i,j,1) = 0.0_rp
210  interp_z2xi_coef( 1:ks-1,i,j,2) = 0.0_rp
211  interp_z2xi_coef( 1:ks-1,i,j,3) = 1.0_rp ! set UNDEF
212 
213  interp_z2xi_idx(ke+1:ka,i,j,1) = ke ! dummmy
214  interp_z2xi_idx(ke+1:ka,i,j,2) = ke ! dummmy
215  interp_z2xi_coef(ke+1:ka,i,j,1) = 0.0_rp
216  interp_z2xi_coef(ke+1:ka,i,j,2) = 0.0_rp
217  interp_z2xi_coef(ke+1:ka,i,j,3) = 1.0_rp ! set UNDEF
218  enddo
219  enddo
220 
221  return
222  end subroutine interp_setup
223 
224  !-----------------------------------------------------------------------------
226  subroutine interp_vertical_xi2z( &
227  var, &
228  var_Z )
229  use scale_const, only: &
231  implicit none
232 
233  real(RP), intent(in) :: var (ka,ia,ja)
234  real(RP), intent(out) :: var_Z(ka,ia,ja)
235 
236  integer :: k, i, j
237  !---------------------------------------------------------------------------
238 
239  do j = jsb, jeb
240  do i = isb, ieb
241  do k = 1, ka
242  var_z(k,i,j) = interp_xi2z_coef(k,i,j,1) * var(interp_xi2z_idx(k,i,j,1),i,j) &
243  + interp_xi2z_coef(k,i,j,2) * var(interp_xi2z_idx(k,i,j,2),i,j) &
244  + interp_xi2z_coef(k,i,j,3) * const_undef
245 
246 ! if ( i == IS .AND. j == JS ) then
247 ! if( IO_L ) write(IO_FID_LOG,*) k,i,j, &
248 ! var_Z(k,i,j), &
249 ! INTERP_xi2z_idx(k,i,j,1), &
250 ! INTERP_xi2z_idx(k,i,j,2), &
251 ! var(INTERP_xi2z_idx(k,i,j,1),i,j), &
252 ! var(INTERP_xi2z_idx(k,i,j,2),i,j)
253 ! endif
254  enddo
255  enddo
256  enddo
257 
258  return
259  end subroutine interp_vertical_xi2z
260 
261  !-----------------------------------------------------------------------------
263  subroutine interp_vertical_z2xi( &
264  var, &
265  var_Xi )
266  use scale_const, only: &
268  implicit none
269 
270  real(RP), intent(in) :: var (ka,ia,ja)
271  real(RP), intent(out) :: var_Xi(ka,ia,ja)
272 
273  integer :: k, i, j
274  !---------------------------------------------------------------------------
275 
276  do j = jsb, jeb
277  do i = isb, ieb
278  do k = 1, ka
279  var_xi(k,i,j) = interp_z2xi_coef(k,i,j,1) * var(interp_z2xi_idx(k,i,j,1),i,j) &
280  + interp_z2xi_coef(k,i,j,2) * var(interp_z2xi_idx(k,i,j,2),i,j) &
281  + interp_z2xi_coef(k,i,j,3) * const_undef
282  enddo
283  enddo
284  enddo
285 
286  return
287  end subroutine interp_vertical_z2xi
288 
289 end module scale_interpolation
290 !-------------------------------------------------------------------------------
subroutine, public interp_setup
Setup.
integer, public jeb
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
subroutine, public interp_vertical_xi2z(var, var_Z)
Reset random seed.
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
real(rp), public const_undef
Definition: scale_const.F90:43
integer, public ieb
module grid index
integer, public ia
of x whole cells (local, with HALO)
module GRID (real space)
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fz
face coordinate [m]: z, local=global
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module GRID (cartesian)
logical, public topo_exist
topography exists?
module profiler
Definition: scale_prof.F90:10
logical, public interp_available
vertical interpolation is available? (have meaning?)
module PRECISION
module TOPOGRAPHY
integer, public isb
subroutine, public interp_vertical_z2xi(var, var_Xi)
Reset random seed.
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public jsb
module INTERPOLATION
integer, public ja
of y whole cells (local, with HALO)