SCALE-RM
Functions/Subroutines | Variables
scale_interpolation Module Reference

module INTERPOLATION More...

Functions/Subroutines

subroutine, public interp_setup
 Setup. More...
 
subroutine, public interp_vertical_xi2z (var, var_Z)
 
subroutine, public interp_vertical_z2xi (var, var_Xi)
 
subroutine, public interp_setup_pres (Kpres)
 Reset random seed. More...
 
subroutine, public interp_update_pres (Kpres, PRES, SFC_PRES, Paxis)
 
subroutine, public interp_vertical_xi2p (Kpres, var, var_P)
 

Variables

logical, public interp_available = .false.
 topography exists & vertical interpolation has meaning? More...
 

Detailed Description

module INTERPOLATION

Description
spacial interpolation module
Author
Team SCALE
History
  • 2013-08-17 (H.Yashiro) [new]

Function/Subroutine Documentation

◆ interp_setup()

subroutine, public scale_interpolation::interp_setup ( )

Setup.

Definition at line 64 of file scale_interpolation.F90.

References scale_grid::grid_cz, scale_grid::grid_fz, scale_grid_index::ia, interp_available, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::ja, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::ks, scale_grid_real::real_cz, scale_grid_real::real_fz, and scale_topography::topo_exist.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

64  use scale_grid, only: &
65  grid_cz, &
66  grid_fz
67  use scale_topography, only: &
69  use scale_grid_real, only: &
70  real_cz, &
71  real_fz
72  implicit none
73 
74  integer :: k, i, j, kk, kp
75  !---------------------------------------------------------------------------
76 
77  if( io_l ) write(io_fid_log,*)
78  if( io_l ) write(io_fid_log,*) '++++++ Module[INTERPOLATION] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
79  if( io_l ) write(io_fid_log,*) '*** No namelists.'
80 
82 
83  if( io_l ) write(io_fid_log,*)
84  if( io_l ) write(io_fid_log,*) '*** Topography exists & interpolation has meaning? : ', interp_available
85 
86  allocate( interp_xi2z_idx(ka,ia,ja,2) )
87  allocate( interp_xi2z_coef(ka,ia,ja,3) )
88  allocate( interp_z2xi_idx(ka,ia,ja,2) )
89  allocate( interp_z2xi_coef(ka,ia,ja,3) )
90 
91  !$omp parallel do default(none) private(i,j,k,kk,kp) OMP_SCHEDULE_ collapse(2) &
92  !$omp shared(JA,IA,KS,KE,GRID_CZ,REAL_FZ,REAL_CZ,INTERP_xi2z_idx,INTERP_xi2z_coef)
93  do j = 1, ja
94  do i = 1, ia
95  do k = ks, ke
96  if ( grid_cz(k) <= real_fz(ks-1,i,j) ) then
97 
98  interp_xi2z_idx(k,i,j,1) = ks ! dummmy
99  interp_xi2z_idx(k,i,j,2) = ks ! dummmy
100  interp_xi2z_coef(k,i,j,1) = 0.0_rp
101  interp_xi2z_coef(k,i,j,2) = 0.0_rp
102  interp_xi2z_coef(k,i,j,3) = 1.0_rp ! set UNDEF
103 
104  elseif( grid_cz(k) <= real_cz(ks,i,j) ) then
105 
106  interp_xi2z_idx(k,i,j,1) = ks ! dummmy
107  interp_xi2z_idx(k,i,j,2) = ks
108  interp_xi2z_coef(k,i,j,1) = 0.0_rp
109  interp_xi2z_coef(k,i,j,2) = 1.0_rp
110  interp_xi2z_coef(k,i,j,3) = 0.0_rp
111 
112  elseif( grid_cz(k) > real_fz(ke,i,j) ) then
113 
114  interp_xi2z_idx(k,i,j,1) = ke ! dummmy
115  interp_xi2z_idx(k,i,j,2) = ke ! dummmy
116  interp_xi2z_coef(k,i,j,1) = 0.0_rp
117  interp_xi2z_coef(k,i,j,2) = 0.0_rp
118  interp_xi2z_coef(k,i,j,3) = 1.0_rp ! set UNDEF
119 
120  elseif( grid_cz(k) > real_cz(ke,i,j) ) then
121 
122  interp_xi2z_idx(k,i,j,1) = ke
123  interp_xi2z_idx(k,i,j,2) = ke ! dummmy
124  interp_xi2z_coef(k,i,j,1) = 1.0_rp
125  interp_xi2z_coef(k,i,j,2) = 0.0_rp
126  interp_xi2z_coef(k,i,j,3) = 0.0_rp
127 
128  else
129 
130  do kk = ks+1, ke
131  kp = kk
132  if( grid_cz(k) <= real_cz(kk,i,j) ) exit
133  enddo
134 
135  interp_xi2z_idx(k,i,j,1) = kp - 1
136  interp_xi2z_idx(k,i,j,2) = kp
137  interp_xi2z_coef(k,i,j,1) = ( real_cz(kp,i,j) - grid_cz(k) ) &
138  / ( real_cz(kp,i,j) - real_cz(kp-1,i,j) )
139  interp_xi2z_coef(k,i,j,2) = ( grid_cz(k) - real_cz(kp-1,i,j) ) &
140  / ( real_cz(kp,i,j) - real_cz(kp-1,i,j) )
141  interp_xi2z_coef(k,i,j,3) = 0.0_rp
142 
143  endif
144  enddo
145  enddo
146  enddo
147 
148  do j = 1, ja
149  do i = 1, ia
150  do k = ks, ke
151  if ( real_cz(k,i,j) <= grid_fz(ks-1) ) then
152 
153  interp_z2xi_idx(k,i,j,1) = ks ! dummmy
154  interp_z2xi_idx(k,i,j,2) = ks ! dummmy
155  interp_z2xi_coef(k,i,j,1) = 0.0_rp
156  interp_z2xi_coef(k,i,j,2) = 0.0_rp
157  interp_z2xi_coef(k,i,j,3) = 1.0_rp ! set UNDEF
158 
159  elseif( real_cz(k,i,j) <= grid_cz(ks) ) then
160 
161  interp_z2xi_idx(k,i,j,1) = ks ! dummmy
162  interp_z2xi_idx(k,i,j,2) = ks
163  interp_z2xi_coef(k,i,j,1) = 0.0_rp
164  interp_z2xi_coef(k,i,j,2) = 1.0_rp
165  interp_z2xi_coef(k,i,j,3) = 0.0_rp
166 
167  elseif( real_cz(k,i,j) > grid_fz(ke) ) then
168 
169  interp_z2xi_idx(k,i,j,1) = ke ! dummmy
170  interp_z2xi_idx(k,i,j,2) = ke ! dummmy
171  interp_z2xi_coef(k,i,j,1) = 0.0_rp
172  interp_z2xi_coef(k,i,j,2) = 0.0_rp
173  interp_z2xi_coef(k,i,j,3) = 1.0_rp ! set UNDEF
174 
175  elseif( real_cz(k,i,j) > grid_cz(ke) ) then
176 
177  interp_z2xi_idx(k,i,j,1) = ke
178  interp_z2xi_idx(k,i,j,2) = ke ! dummmy
179  interp_z2xi_coef(k,i,j,1) = 1.0_rp
180  interp_z2xi_coef(k,i,j,2) = 0.0_rp
181  interp_z2xi_coef(k,i,j,3) = 0.0_rp
182 
183  else
184 
185  do kk = ks+1, ke
186  kp = kk
187  if( real_cz(k,i,j) <= grid_cz(kk) ) exit
188  enddo
189 
190  interp_z2xi_idx(k,i,j,1) = kp - 1
191  interp_z2xi_idx(k,i,j,2) = kp
192  interp_z2xi_coef(k,i,j,1) = ( grid_cz(kp) - real_cz(k,i,j) ) &
193  / ( grid_cz(kp) - grid_cz(kp-1) )
194  interp_z2xi_coef(k,i,j,2) = ( real_cz(k,i,j) - grid_cz(kp-1) ) &
195  / ( grid_cz(kp) - grid_cz(kp-1) )
196  interp_z2xi_coef(k,i,j,3) = 0.0_rp
197 
198  endif
199  enddo
200  enddo
201  enddo
202 
203  !$omp parallel do default(none) private(i,j) OMP_SCHEDULE_ collapse(2) &
204  !$omp shared(JA,IA,INTERP_xi2z_idx,INTERP_xi2z_coef,INTERP_z2xi_idx,INTERP_z2xi_coef,KS,KE,KA)
205  do j = 1, ja
206  do i = 1, ia
207  interp_xi2z_idx( 1:ks-1,i,j,1) = ks ! dummmy
208  interp_xi2z_idx( 1:ks-1,i,j,2) = ks ! dummmy
209  interp_xi2z_coef( 1:ks-1,i,j,1) = 0.0_rp
210  interp_xi2z_coef( 1:ks-1,i,j,2) = 0.0_rp
211  interp_xi2z_coef( 1:ks-1,i,j,3) = 1.0_rp ! set UNDEF
212 
213  interp_xi2z_idx(ke+1:ka,i,j,1) = ke ! dummmy
214  interp_xi2z_idx(ke+1:ka,i,j,2) = ke ! dummmy
215  interp_xi2z_coef(ke+1:ka,i,j,1) = 0.0_rp
216  interp_xi2z_coef(ke+1:ka,i,j,2) = 0.0_rp
217  interp_xi2z_coef(ke+1:ka,i,j,3) = 1.0_rp ! set UNDEF
218 
219  interp_z2xi_idx( 1:ks-1,i,j,1) = ks ! dummmy
220  interp_z2xi_idx( 1:ks-1,i,j,2) = ks ! dummmy
221  interp_z2xi_coef( 1:ks-1,i,j,1) = 0.0_rp
222  interp_z2xi_coef( 1:ks-1,i,j,2) = 0.0_rp
223  interp_z2xi_coef( 1:ks-1,i,j,3) = 1.0_rp ! set UNDEF
224 
225  interp_z2xi_idx(ke+1:ka,i,j,1) = ke ! dummmy
226  interp_z2xi_idx(ke+1:ka,i,j,2) = ke ! dummmy
227  interp_z2xi_coef(ke+1:ka,i,j,1) = 0.0_rp
228  interp_z2xi_coef(ke+1:ka,i,j,2) = 0.0_rp
229  interp_z2xi_coef(ke+1:ka,i,j,3) = 1.0_rp ! set UNDEF
230  enddo
231  enddo
232 
233  return
module GRID (real space)
module GRID (cartesian)
logical, public topo_exist
topography exists?
logical, public interp_available
topography exists & vertical interpolation has meaning?
module TOPOGRAPHY
Here is the caller graph for this function:

◆ interp_vertical_xi2z()

subroutine, public scale_interpolation::interp_vertical_xi2z ( real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Z 
)

Definition at line 240 of file scale_interpolation.F90.

References scale_const::const_undef, scale_grid::grid_cz, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, and scale_grid_real::real_cz.

Referenced by scale_atmos_refstate::atmos_refstate_update(), and scale_history::hist_put_3d().

240  use scale_const, only: &
241  const_undef
242  use scale_grid, only: &
243  grid_cz
244  use scale_grid_real, only: &
245  real_cz
246  use scale_matrix, only: &
247  matrix_solver_tridiagonal
248  implicit none
249 
250  real(RP), intent(in) :: var (KA,IA,JA)
251  real(RP), intent(out) :: var_Z(KA,IA,JA)
252 
253  real(RP) :: FDZ(KA)
254  real(RP) :: MD(KA)
255  real(RP) :: U(KA)
256  real(RP) :: V(KA)
257  real(RP) :: c1, c2, c3
258  real(RP) :: d
259  integer :: kk
260 
261  integer :: k, i, j
262  !---------------------------------------------------------------------------
263 
264  if ( kmax == 2 ) then
265 
266  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
267  !$omp private(k,i,j) &
268  !$omp shared(KA,ISB,IEB,JSB,JEB) &
269  !$omp shared(var,var_Z) &
270  !$omp shared(INTERP_xi2z_idx,INTERP_xi2z_coef,CONST_UNDEF)
271  do j = jsb, jeb
272  do i = isb, ieb
273  do k = 1, ka
274  var_z(k,i,j) = interp_xi2z_coef(k,i,j,1) * var(interp_xi2z_idx(k,i,j,1),i,j) &
275  + interp_xi2z_coef(k,i,j,2) * var(interp_xi2z_idx(k,i,j,2),i,j) &
276  + interp_xi2z_coef(k,i,j,3) * const_undef
277  enddo
278  enddo
279  enddo
280 
281  else
282 
283  ! cubic spline
284  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
285  !$omp private(k,i,j) &
286  !$omp private(FDZ,MD,U,V,kk,c1,c2,c3,d) &
287  !$omp shared(KA,KS,KE,KMAX,ISB,IEB,JSB,JEB) &
288  !$omp shared(var,var_Z) &
289  !$omp shared(INTERP_xi2z_idx,INTERP_xi2z_coef,GRID_CZ,REAL_CZ,CONST_UNDEF)
290  do j = jsb, jeb
291  do i = isb, ieb
292 
293  do k = ks, ke-1
294  fdz(k) = real_cz(k+1,i,j) - real_cz(k,i,j)
295  end do
296 
297  md(ks+1) = 2.0 * ( fdz(ks) + fdz(ks+1) ) + fdz(ks)
298  do k = ks+2, ke-2
299  md(k) = 2.0 * ( fdz(k-1) + fdz(k) )
300  end do
301  md(ke-1) = 2.0 * ( fdz(ke-2) + fdz(ke-1) ) + fdz(ke-1)
302 
303  do k = ks+1, ke-1
304  v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / fdz(k) &
305  - ( var(k ,i,j) - var(k-1,i,j) ) / fdz(k-1)
306  end do
307 
308  call matrix_solver_tridiagonal( &
309  kmax-2, &
310  fdz(ks+1:ke-1), md(ks+1:ke-1), fdz(ks+1:ke-1), v(ks+1:ke-1), & ! (in)
311  u(ks+1:ke-1) ) ! (out)
312  u(ks) = u(ks+1)
313  u(ke) = u(ke-1)
314 
315  do k = 1, ka
316  kk = min(interp_xi2z_idx(k,i,j,1), ke-1)
317  c3 = ( u(kk+1) - u(kk) ) / fdz(kk)
318  c2 = 3.0_rp * u(kk)
319  c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / fdz(kk) - ( u(kk) * 2.0_rp + u(kk+1) ) * fdz(kk)
320  d = grid_cz(k) - real_cz(kk,i,j)
321  var_z(k,i,j) = interp_xi2z_coef(k,i,j,3) * const_undef &
322  + ( 1.0_rp - interp_xi2z_coef(k,i,j,3) ) &
323  * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
324  end do
325 
326  end do
327  end do
328 
329  end if
330 
331  return
module GRID (real space)
module MATRIX
module CONSTANT
Definition: scale_const.F90:14
module GRID (cartesian)
Here is the caller graph for this function:

◆ interp_vertical_z2xi()

subroutine, public scale_interpolation::interp_vertical_z2xi ( real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Xi 
)

Definition at line 338 of file scale_interpolation.F90.

References scale_const::const_undef, scale_grid::grid_cz, scale_grid::grid_fdz, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, and scale_grid_real::real_cz.

Referenced by scale_atmos_refstate::atmos_refstate_calc3d().

338  use scale_const, only: &
339  const_undef
340  use scale_grid, only: &
341  grid_cz, &
342  fdz => grid_fdz
343  use scale_matrix, only: &
344  matrix_solver_tridiagonal
345  use scale_grid_real, only: &
346  real_cz
347  implicit none
348 
349  real(RP), intent(in) :: var (KA,IA,JA)
350  real(RP), intent(out) :: var_Xi(KA,IA,JA)
351 
352  real(RP) :: MD(KA)
353  real(RP) :: U(KA)
354  real(RP) :: V(KA)
355  real(RP) :: c1, c2, c3
356  real(RP) :: d
357  integer :: kk
358 
359  integer :: k, i, j
360  !---------------------------------------------------------------------------
361 
362  if ( kmax == 2 ) then
363 
364  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
365  !$omp private(k,i,j) &
366  !$omp shared(KA,ISB,IEB,JSB,JEB) &
367  !$omp shared(var,var_Xi) &
368  !$omp shared(INTERP_z2xi_idx,INTERP_z2xi_coef,CONST_UNDEF)
369  do j = jsb, jeb
370  do i = isb, ieb
371  do k = 1, ka
372  var_xi(k,i,j) = interp_z2xi_coef(k,i,j,1) * var(interp_z2xi_idx(k,i,j,1),i,j) &
373  + interp_z2xi_coef(k,i,j,2) * var(interp_z2xi_idx(k,i,j,2),i,j) &
374  + interp_z2xi_coef(k,i,j,3) * const_undef
375  enddo
376  enddo
377  enddo
378 
379  else
380 
381  ! cubic spline
382 
383  md(ks+1) = 2.0 * ( fdz(ks) + fdz(ks+1) ) + fdz(ks)
384  do k = ks+2, ke-2
385  md(k) = 2.0 * ( fdz(k-1) + fdz(k) )
386  end do
387  md(ke-1) = 2.0 * ( fdz(ke-2) + fdz(ke-1) ) + fdz(ke-1)
388 
389  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
390  !$omp private(k,i,j) &
391  !$omp private(U,V,kk,c1,c2,c3,d) &
392  !$omp shared(KA,KS,KE,KMAX,ISB,IEB,JSB,JEB) &
393  !$omp shared(var,var_Xi,MD,FDZ) &
394  !$omp shared(INTERP_z2xi_idx,INTERP_z2xi_coef,GRID_CZ,REAL_CZ,CONST_UNDEF)
395  do j = jsb, jeb
396  do i = isb, ieb
397 
398  do k = ks+1, ke-1
399  v(k) = ( var(k+1,i,j) - var(k ,i,j) ) / fdz(k) &
400  - ( var(k ,i,j) - var(k-1,i,j) ) / fdz(k-1)
401  end do
402 
403  call matrix_solver_tridiagonal( &
404  kmax-2, &
405  fdz(ks+1:ke-1), md(ks+1:ke-1), fdz(ks+1:ke-1), v(ks+1:ke-1), & ! (in)
406  u(ks+1:ke-1) ) ! (out)
407  u(ks) = u(ks+1)
408  u(ke) = u(ke-1)
409 
410  do k = 1, ka
411  kk = min(interp_z2xi_idx(k,i,j,1), ke-1)
412  c3 = ( u(kk+1) - u(kk) ) / fdz(kk)
413  c2 = 3.0_rp * u(kk)
414  c1 = ( var(kk+1,i,j) - var(kk,i,j) ) / fdz(kk) - ( u(kk) * 2.0_rp + u(kk+1) ) * fdz(kk)
415  d = real_cz(k,i,j) - grid_cz(kk)
416  var_xi(k,i,j) = interp_z2xi_coef(k,i,j,3) * const_undef &
417  + ( 1.0_rp - interp_z2xi_coef(k,i,j,3) ) &
418  * ( ( ( c3 * d + c2 ) * d + c1 ) * d + var(kk,i,j) )
419  end do
420 
421  end do
422  end do
423 
424  end if
425 
426  return
module GRID (real space)
module MATRIX
module CONSTANT
Definition: scale_const.F90:14
module GRID (cartesian)
Here is the caller graph for this function:

◆ interp_setup_pres()

subroutine, public scale_interpolation::interp_setup_pres ( integer, intent(in)  Kpres)

Reset random seed.

Definition at line 433 of file scale_interpolation.F90.

References scale_grid_index::ia, and scale_grid_index::ja.

Referenced by scale_history::hist_setup().

433  implicit none
434 
435  integer, intent(in) :: Kpres
436  !---------------------------------------------------------------------------
437 
438  allocate( interp_xi2p_idx(kpres,ia,ja,2) )
439  allocate( interp_xi2p_coef(kpres,ia,ja,3) )
440 
441  return
Here is the caller graph for this function:

◆ interp_update_pres()

subroutine, public scale_interpolation::interp_update_pres ( integer, intent(in)  Kpres,
real(rp), dimension (ka,ia,ja), intent(in)  PRES,
real(rp), dimension( ia,ja), intent(in)  SFC_PRES,
real(rp), dimension (kpres), intent(in)  Paxis 
)

Definition at line 450 of file scale_interpolation.F90.

References scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, scale_grid_index::jsb, scale_grid_index::ke, and scale_grid_index::ks.

Referenced by scale_history::hist_setpres().

450  implicit none
451 
452  integer, intent(in) :: Kpres
453  real(RP), intent(in) :: PRES (KA,IA,JA) ! pressure in Xi coordinate [Pa]
454  real(RP), intent(in) :: SFC_PRES( IA,JA) ! surface pressure [Pa]
455  real(RP), intent(in) :: Paxis (Kpres) ! pressure level to output [Pa]
456 
457  real(RP) :: LnPRES (KA,IA,JA) ! (log) pressure in Xi coordinate [Pa]
458  real(RP) :: LnSFC_PRES( IA,JA) ! (log) surface pressure [Pa]
459  real(RP) :: LnPaxis (Kpres) ! (log) pressure level to output [Pa]
460 
461  integer :: k, i, j, kk, kp
462  !---------------------------------------------------------------------------
463 
464  do j = jsb, jeb
465  do i = isb, ieb
466  do k = ks, ke
467  lnpres(k,i,j) = log( pres(k,i,j) )
468  enddo
469  enddo
470  enddo
471 
472  do j = jsb, jeb
473  do i = isb, ieb
474  lnsfc_pres(i,j) = log( sfc_pres(i,j) )
475  enddo
476  enddo
477 
478  lnpaxis(:) = log( paxis(:) )
479 
480  do j = jsb, jeb
481  do i = isb, ieb
482  do k = 1, kpres
483  if ( lnpaxis(k) >= lnsfc_pres(i,j) ) then
484 
485  interp_xi2p_idx(k,i,j,1) = ks ! dummmy
486  interp_xi2p_idx(k,i,j,2) = ks ! dummmy
487  interp_xi2p_coef(k,i,j,1) = 0.0_rp
488  interp_xi2p_coef(k,i,j,2) = 0.0_rp
489  interp_xi2p_coef(k,i,j,3) = 1.0_rp ! set UNDEF
490 
491  elseif( lnpaxis(k) >= lnpres(ks,i,j) ) then
492 
493  interp_xi2p_idx(k,i,j,1) = ks ! dummmy
494  interp_xi2p_idx(k,i,j,2) = ks
495  interp_xi2p_coef(k,i,j,1) = 0.0_rp
496  interp_xi2p_coef(k,i,j,2) = 1.0_rp
497  interp_xi2p_coef(k,i,j,3) = 0.0_rp
498 
499  elseif( lnpaxis(k) < lnpres(ke,i,j) ) then
500 
501  interp_xi2p_idx(k,i,j,1) = ke ! dummmy
502  interp_xi2p_idx(k,i,j,2) = ke ! dummmy
503  interp_xi2p_coef(k,i,j,1) = 0.0_rp
504  interp_xi2p_coef(k,i,j,2) = 0.0_rp
505  interp_xi2p_coef(k,i,j,3) = 1.0_rp ! set UNDEF
506 
507  else
508 
509  do kk = ks+1, ke
510  kp = kk
511  if( lnpaxis(k) >= lnpres(kk,i,j) ) exit
512  enddo
513 
514  interp_xi2p_idx(k,i,j,1) = kp - 1
515  interp_xi2p_idx(k,i,j,2) = kp
516  interp_xi2p_coef(k,i,j,1) = ( lnpaxis(k) - lnpres(kp,i,j) ) &
517  / ( lnpres(kp-1,i,j) - lnpres(kp,i,j) )
518  interp_xi2p_coef(k,i,j,2) = ( lnpres(kp-1,i,j) - lnpaxis(k) ) &
519  / ( lnpres(kp-1,i,j) - lnpres(kp,i,j) )
520  interp_xi2p_coef(k,i,j,3) = 0.0_rp
521 
522  endif
523  enddo
524  enddo
525  enddo
526 
527  return
Here is the caller graph for this function:

◆ interp_vertical_xi2p()

subroutine, public scale_interpolation::interp_vertical_xi2p ( integer, intent(in)  Kpres,
real(rp), dimension (ka ,ia,ja), intent(in)  var,
real(rp), dimension(kpres,ia,ja), intent(out)  var_P 
)

Definition at line 535 of file scale_interpolation.F90.

References scale_const::const_undef, scale_grid_index::ieb, scale_grid_index::isb, scale_grid_index::jeb, and scale_grid_index::jsb.

Referenced by scale_history::hist_put_3d().

535  use scale_const, only: &
536  const_undef
537  implicit none
538 
539  integer, intent(in) :: Kpres
540  real(RP), intent(in) :: var (KA ,IA,JA)
541  real(RP), intent(out) :: var_P(Kpres,IA,JA)
542 
543  integer :: k, i, j
544  !---------------------------------------------------------------------------
545 
546  do j = jsb, jeb
547  do i = isb, ieb
548  do k = 1, kpres
549  var_p(k,i,j) = interp_xi2p_coef(k,i,j,1) * var(interp_xi2p_idx(k,i,j,1),i,j) &
550  + interp_xi2p_coef(k,i,j,2) * var(interp_xi2p_idx(k,i,j,2),i,j) &
551  + interp_xi2p_coef(k,i,j,3) * const_undef
552  enddo
553  enddo
554  enddo
555 
556  return
module CONSTANT
Definition: scale_const.F90:14
Here is the caller graph for this function:

Variable Documentation

◆ interp_available

logical, public scale_interpolation::interp_available = .false.

topography exists & vertical interpolation has meaning?

Definition at line 41 of file scale_interpolation.F90.

Referenced by scale_history::hist_put_3d(), and interp_setup().

41  logical, public :: INTERP_available = .false.