SCALE-RM
Functions/Subroutines | Variables
scale_interp_vert Module Reference

module INTERPOLATION vertical More...

Functions/Subroutines

subroutine, public interp_vert_setcoef (KA, KS, KE, IA, IS, IE, JA, JS, JE, TOPO_exist, Xi, Xih, Z, Zh)
 Setup. More...
 
subroutine, public interp_vert_xi2z (KA, KS, KE, IA, IS, IE, JA, JS, JE, Xi, Z, var, var_Z)
 
subroutine, public interp_vert_z2xi (KA, KS, KE, IA, IS, IE, JA, JS, JE, Z, Xi, var, var_Xi)
 
subroutine, public interp_vert_xih2zh (KA, KS, KE, IA, IS, IE, JA, JS, JE, Xih, Zh, var, var_Z)
 
subroutine, public interp_vert_zh2xih (KA, KS, KE, IA, IS, IE, JA, JS, JE, Zh, Xih, var, var_Xi)
 
subroutine, public interp_vert_alloc_pres (Kpres, KA, IA, JA)
 Setup. More...
 
subroutine, public interp_vert_setcoef_pres (Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, PRES, PRESh, SFC_PRES, Paxis)
 
subroutine, public interp_vert_xi2p (Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
 
subroutine, public interp_vert_xih2p (Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
 

Variables

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

Detailed Description

module INTERPOLATION vertical

Description
spacial interpolation module, vertical only (xi<->z,xi->pres)
Author
Team SCALE

Function/Subroutine Documentation

◆ interp_vert_setcoef()

subroutine, public scale_interp_vert::interp_vert_setcoef ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
logical, intent(in)  TOPO_exist,
real(rp), dimension ( ka), intent(in)  Xi,
real(rp), dimension(0:ka), intent(in)  Xih,
real(rp), dimension ( ka,ia,ja), intent(in)  Z,
real(rp), dimension (0:ka,ia,ja), intent(in)  Zh 
)

Setup.

Definition at line 86 of file scale_interp_vert.F90.

86  use scale_interp, only: &
88  implicit none
89  integer, intent(in) :: KA, KS, KE
90  integer, intent(in) :: IA, IS, IE
91  integer, intent(in) :: JA, JS, JE
92  logical, intent(in) :: TOPO_exist
93  real(RP), intent(in) :: Xi ( KA)
94  real(RP), intent(in) :: Xih(0:KA)
95  real(RP), intent(in) :: Z ( KA,IA,JA)
96  real(RP), intent(in) :: Zh (0:KA,IA,JA)
97 
98  integer :: i, j
99  !---------------------------------------------------------------------------
100 
101  log_newline
102  log_info("INTERP_VERT_setcoef",*) 'Setup'
103  log_info("INTERP_VERT_setcoef",*) 'No namelists.'
104 
105  interp_available = topo_exist
106 
107  log_newline
108  log_info("INTERP_VERT_setcoef",*) 'Topography exists & interpolation has meaning? : ', interp_available
109 
110  ! full level
111 
112  allocate( interp_xi2z_idx(ka,2,ia,ja) )
113  allocate( interp_xi2z_coef(ka, ia,ja) )
114  allocate( interp_z2xi_idx(ka,2,ia,ja) )
115  allocate( interp_z2xi_coef(ka, ia,ja) )
116 
117  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
118  !$omp shared(JA,IA,KA,KS,KE,Xi,Z,INTERP_xi2z_idx,INTERP_xi2z_coef)
119  do j = 1, ja
120  do i = 1, ia
121  call interp_factor1d( ka, ks, ke, ka, ks, ke, &
122  z(:,i,j), xi(:), & ! (in)
123  interp_xi2z_idx(:,:,i,j), & ! (in)
124  interp_xi2z_coef(:, i,j), & ! (out)
125  flag_extrap = .true. ) ! (in)
126  enddo
127  enddo
128 
129  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
130  !$omp shared(JA,IA,KA,KS,KE,Z,Xi,INTERP_z2xi_idx,INTERP_z2xi_coef)
131  do j = 1, ja
132  do i = 1, ia
133  call interp_factor1d( ka, ks, ke, ka, ks, ke, &
134  xi(:), z(:,i,j), & ! (in)
135  interp_z2xi_idx(:,:,i,j), & ! (in)
136  interp_z2xi_coef(:, i,j), & ! (out)
137  flag_extrap = .true. ) ! (in)
138  enddo
139  enddo
140 
141 
142  ! half level
143 
144  allocate( interp_xih2zh_idx(ka,2,ia,ja) )
145  allocate( interp_xih2zh_coef(ka, ia,ja) )
146  allocate( interp_zh2xih_idx(ka,2,ia,ja) )
147  allocate( interp_zh2xih_coef(ka, ia,ja) )
148 
149  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
150  !$omp shared(JA,IA,KA,KS,KE,Xih,Zh,INTERP_xih2zh_idx,INTERP_xih2zh_coef)
151  do j = 1, ja
152  do i = 1, ia
153  call interp_factor1d( ka, ks-1, ke, ka, ks-1, ke, &
154  zh(1:,i,j), xih(1:), & ! (in)
155  interp_xih2zh_idx(:,:,i,j), & ! (in)
156  interp_xih2zh_coef(:, i,j), & ! (out)
157  flag_extrap = .false. ) ! (in)
158  enddo
159  enddo
160 
161  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
162  !$omp shared(JA,IA,KA,KS,KE,Xih,Zh,INTERP_zh2xih_idx,INTERP_zh2xih_coef)
163  do j = 1, ja
164  do i = 1, ia
165  call interp_factor1d( ka, ks-1, ke, ka, ks-1, ke, &
166  xih(1:), zh(1:,i,j), & ! (in)
167  interp_zh2xih_idx(:,:,i,j), & ! (in)
168  interp_zh2xih_coef(:, i,j), & ! (out)
169  flag_extrap = .true. ) ! (in)
170  enddo
171  enddo
172 
173 
174  return

References interp_available, and scale_interp::interp_factor1d().

Referenced by scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_setup().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_xi2z()

subroutine, public scale_interp_vert::interp_vert_xi2z ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka), intent(in)  Xi,
real(rp), dimension (ka,ia,ja), intent(in)  Z,
real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Z 
)

Definition at line 186 of file scale_interp_vert.F90.

186  use scale_interp, only: &
188  implicit none
189  integer, intent(in) :: KA, KS, KE
190  integer, intent(in) :: IA, IS, IE
191  integer, intent(in) :: JA, JS, JE
192  real(RP), intent(in) :: Xi (KA)
193  real(RP), intent(in) :: Z (KA,IA,JA)
194  real(RP), intent(in) :: var (KA,IA,JA)
195  real(RP), intent(out) :: var_Z(KA,IA,JA)
196 
197  integer :: i, j
198  !---------------------------------------------------------------------------
199 
200  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
201  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
202  !$omp shared(Xi,Z,var,var_Z) &
203  !$omp shared(INTERP_xi2z_idx,INTERP_xi2z_coef)
204  do j = js, je
205  do i = is, ie
206  call interp_interp1d( ka, ks, ke, ka, ks, ke, &
207  interp_xi2z_idx(:,:,i,j), & ! (in)
208  interp_xi2z_coef(:, i,j), & ! (in)
209  z(:,i,j), xi(:), & ! (in)
210  var(:,i,j), & ! (in)
211  var_z(:,i,j) ) ! (out)
212  end do
213  end do
214 
215  return

References scale_interp::interp_interp1d().

Referenced by scale_atmos_refstate::atmos_refstate_update(), and scale_file_history_cartesc::file_history_cartesc_truncate_3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_z2xi()

subroutine, public scale_interp_vert::interp_vert_z2xi ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka,ia,ja), intent(in)  Z,
real(rp), dimension (ka), intent(in)  Xi,
real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Xi 
)

Definition at line 227 of file scale_interp_vert.F90.

227  use scale_interp, only: &
229  implicit none
230  integer, intent(in) :: KA, KS, KE
231  integer, intent(in) :: IA, IS, IE
232  integer, intent(in) :: JA, JS, JE
233  real(RP), intent(in) :: Z (KA,IA,JA)
234  real(RP), intent(in) :: Xi (KA)
235  real(RP), intent(in) :: var (KA,IA,JA)
236  real(RP), intent(out) :: var_Xi(KA,IA,JA)
237 
238  integer :: i, j
239  !---------------------------------------------------------------------------
240 
241  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
242  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
243  !$omp shared(Z,Xi,var,var_Xi) &
244  !$omp shared(INTERP_z2xi_idx,INTERP_z2xi_coef)
245  do j = js, je
246  do i = is, ie
247  call interp_interp1d( ka, ks, ke, ka, ks, ke, &
248  interp_z2xi_idx(:,:,i,j), & ! (in)
249  interp_z2xi_coef(:, i,j), & ! (in)
250  xi(:), z(:,i,j), & ! (in)
251  var(:,i,j), & ! (in)
252  var_xi(:,i,j) ) ! (out)
253  end do
254  end do
255 
256  return

References scale_interp::interp_interp1d().

Referenced by scale_atmos_refstate::atmos_refstate_calc3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_xih2zh()

subroutine, public scale_interp_vert::interp_vert_xih2zh ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (0:ka), intent(in)  Xih,
real(rp), dimension (0:ka,ia,ja), intent(in)  Zh,
real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Z 
)

Definition at line 268 of file scale_interp_vert.F90.

268  use scale_interp, only: &
270  implicit none
271  integer, intent(in) :: KA, KS, KE
272  integer, intent(in) :: IA, IS, IE
273  integer, intent(in) :: JA, JS, JE
274  real(RP), intent(in) :: Xih (0:KA)
275  real(RP), intent(in) :: Zh (0:KA,IA,JA)
276  real(RP), intent(in) :: var (KA,IA,JA)
277  real(RP), intent(out) :: var_Z(KA,IA,JA)
278 
279  integer :: i, j
280  !---------------------------------------------------------------------------
281 
282  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
283  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
284  !$omp shared(Xih,Zh,var,var_Z) &
285  !$omp shared(INTERP_xih2zh_idx,INTERP_xih2zh_coef)
286  do j = js, je
287  do i = is, ie
288  call interp_interp1d( ka, ks-1, ke, ka, ks-1, ke, &
289  interp_xih2zh_idx(:,:,i,j), & ! (in)
290  interp_xih2zh_coef(:, i,j), & ! (in)
291  zh(1:,i,j), xih(1:), & ! (in)
292  var(:,i,j), & ! (in)
293  var_z(:,i,j) ) ! (out)
294  enddo
295  enddo
296 
297  return

References scale_interp::interp_interp1d().

Referenced by scale_file_history_cartesc::file_history_cartesc_truncate_3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_zh2xih()

subroutine, public scale_interp_vert::interp_vert_zh2xih ( integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (0:ka,ia,ja), intent(in)  Zh,
real(rp), dimension (0:ka), intent(in)  Xih,
real(rp), dimension (ka,ia,ja), intent(in)  var,
real(rp), dimension(ka,ia,ja), intent(out)  var_Xi 
)

Definition at line 309 of file scale_interp_vert.F90.

309  use scale_interp, only: &
311  implicit none
312  integer, intent(in) :: KA, KS, KE
313  integer, intent(in) :: IA, IS, IE
314  integer, intent(in) :: JA, JS, JE
315  real(RP), intent(in) :: Zh (0:KA,IA,JA)
316  real(RP), intent(in) :: Xih (0:KA)
317  real(RP), intent(in) :: var (KA,IA,JA)
318  real(RP), intent(out) :: var_Xi(KA,IA,JA)
319 
320  integer :: i, j
321  !---------------------------------------------------------------------------
322 
323  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
324  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
325  !$omp shared(Zh,Xih,var,var_Xi) &
326  !$omp shared(INTERP_zh2xih_idx,INTERP_zh2xih_coef)
327  do j = js, je
328  do i = is, ie
329  call interp_interp1d( ka, ks-1, ke, ka, ks-1, ke, &
330  interp_zh2xih_idx(:,:,i,j), & ! (in)
331  interp_zh2xih_coef(:, i,j), & ! (in)
332  xih(1:), zh(1:,i,j), & ! (in)
333  var(:,i,j), & ! (in)
334  var_xi(:,i,j) ) ! (out)
335  enddo
336  enddo
337 
338  return

References scale_interp::interp_interp1d().

Here is the call graph for this function:

◆ interp_vert_alloc_pres()

subroutine, public scale_interp_vert::interp_vert_alloc_pres ( integer, intent(in)  Kpres,
integer, intent(in)  KA,
integer, intent(in)  IA,
integer, intent(in)  JA 
)

Setup.

Definition at line 345 of file scale_interp_vert.F90.

345  implicit none
346  integer, intent(in) :: Kpres
347  integer, intent(in) :: KA
348  integer, intent(in) :: IA
349  integer, intent(in) :: JA
350  !---------------------------------------------------------------------------
351 
352  allocate( interp_xi2p_idx(kpres,2,ia,ja) )
353  allocate( interp_xi2p_coef(kpres, ia,ja) )
354  allocate( interp_xih2p_idx(kpres,2,ia,ja) )
355  allocate( interp_xih2p_coef(kpres, ia,ja) )
356 
357  allocate( lnpres(ka,ia,ja) )
358  allocate( lnpresh(ka,ia,ja) )
359  allocate( lnpaxis(kpres) )
360 
361  return

Referenced by scale_file_history_cartesc::file_history_cartesc_setup().

Here is the caller graph for this function:

◆ interp_vert_setcoef_pres()

subroutine, public scale_interp_vert::interp_vert_setcoef_pres ( integer, intent(in)  Kpres,
integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension ( ka,ia,ja), intent(in)  PRES,
real(rp), dimension (0:ka,ia,ja), intent(in)  PRESh,
real(rp), dimension( ia,ja), intent(in)  SFC_PRES,
real(rp), dimension (kpres), intent(in)  Paxis 
)

Definition at line 374 of file scale_interp_vert.F90.

374  use scale_interp, only: &
376  implicit none
377 
378  integer, intent(in) :: Kpres
379  integer, intent(in) :: KA, KS, KE
380  integer, intent(in) :: IA, IS, IE
381  integer, intent(in) :: JA, JS, JE
382  real(RP), intent(in) :: PRES ( KA,IA,JA) ! pressure in Xi coordinate [Pa]
383  real(RP), intent(in) :: PRESh (0:KA,IA,JA) ! pressure in Xi coordinate [Pa], layer interface
384  real(RP), intent(in) :: SFC_PRES( IA,JA) ! surface pressure [Pa]
385  real(RP), intent(in) :: Paxis (Kpres) ! pressure level to output [Pa]
386 
387  integer :: k, i, j
388  !---------------------------------------------------------------------------
389 
390  ! full level
391 
392 !OCL SERIAL
393  do k = 1, kpres
394  lnpaxis(k) = - log( paxis(k) )
395  end do
396 
397  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
398  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
399  !$omp shared(PRES,LnPRES,LnPaxis) &
400  !$omp shared(INTERP_xi2p_idx,INTERP_xi2p_coef)
401  do j = js, je
402  do i = is, ie
403  do k = ks, ke
404  lnpres(k,i,j) = - log( pres(k,i,j) )
405  end do
406  call interp_factor1d( ka, ks, ke, kpres, 1, kpres, &
407  lnpres(:,i,j), lnpaxis(:), & ! (in)
408  interp_xi2p_idx(:,:,i,j), & ! (in)
409  interp_xi2p_coef(:, i,j), & ! (out)
410  flag_extrap = .false. ) ! (in)
411  enddo
412  enddo
413 
414  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
415  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
416  !$omp shared(PRESh,LnPRESh,LnPaxis) &
417  !$omp shared(INTERP_xih2p_idx,INTERP_xih2p_coef)
418  do j = js, je
419  do i = is, ie
420  do k = ks-1, ke
421  lnpresh(k,i,j) = - log( presh(k,i,j) )
422  end do
423  call interp_factor1d( ka, ks-1, ke, kpres, 1, kpres, &
424  lnpresh(:,i,j), lnpaxis(:), & ! (in)
425  interp_xih2p_idx(:,:,i,j), & ! (in)
426  interp_xih2p_coef(:, i,j), & ! (out)
427  flag_extrap = .false. ) ! (in)
428  enddo
429  enddo
430 
431 
432  return

References scale_interp::interp_factor1d().

Referenced by scale_file_history_cartesc::file_history_cartesc_set_pres().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_xi2p()

subroutine, public scale_interp_vert::interp_vert_xi2p ( integer, intent(in)  Kpres,
integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka ,ia,ja), intent(in)  var,
real(rp), dimension(kpres,ia,ja), intent(out)  var_P 
)

Definition at line 443 of file scale_interp_vert.F90.

443  use scale_interp, only: &
445  implicit none
446  integer, intent(in) :: Kpres
447  integer, intent(in) :: KA, KS, KE
448  integer, intent(in) :: IA, IS, IE
449  integer, intent(in) :: JA, JS, JE
450  real(RP), intent(in) :: var (KA ,IA,JA)
451  real(RP), intent(out) :: var_P(Kpres,IA,JA)
452 
453  integer :: i, j
454  !---------------------------------------------------------------------------
455 
456  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
457  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
458  !$omp shared(LnPRES,LnPaxis,var,var_P) &
459  !$omp shared(INTERP_xi2p_idx,INTERP_xi2p_coef)
460  do j = js, je
461  do i = is, ie
462  call interp_interp1d( ka, ks, ke, kpres, 1, kpres, &
463  interp_xi2p_idx(:,:,i,j), & ! (in)
464  interp_xi2p_coef(:, i,j), & ! (in)
465  lnpres(:,i,j), lnpaxis(:), & ! (in)
466  var(:,i,j), & ! (in)
467  var_p(:,i,j) ) ! (out)
468  enddo
469  enddo
470 
471  return

References scale_interp::interp_interp1d().

Referenced by scale_file_history_cartesc::file_history_cartesc_truncate_3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ interp_vert_xih2p()

subroutine, public scale_interp_vert::interp_vert_xih2p ( integer, intent(in)  Kpres,
integer, intent(in)  KA,
integer, intent(in)  KS,
integer, intent(in)  KE,
integer, intent(in)  IA,
integer, intent(in)  IS,
integer, intent(in)  IE,
integer, intent(in)  JA,
integer, intent(in)  JS,
integer, intent(in)  JE,
real(rp), dimension (ka ,ia,ja), intent(in)  var,
real(rp), dimension(kpres,ia,ja), intent(out)  var_P 
)

Definition at line 482 of file scale_interp_vert.F90.

482  use scale_interp, only: &
484  implicit none
485 
486  integer, intent(in) :: Kpres
487  integer, intent(in) :: KA, KS, KE
488  integer, intent(in) :: IA, IS, IE
489  integer, intent(in) :: JA, JS, JE
490  real(RP), intent(in) :: var (KA ,IA,JA)
491  real(RP), intent(out) :: var_P(Kpres,IA,JA)
492 
493  integer :: k, i, j
494  !---------------------------------------------------------------------------
495 
496  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
497  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
498  !$omp shared(LnPRESh,LnPaxis,var,var_P) &
499  !$omp shared(INTERP_xih2p_idx,INTERP_xih2p_coef)
500  do j = js, je
501  do i = is, ie
502  call interp_interp1d( ka, ks-1, ke, kpres, 1, kpres, &
503  interp_xih2p_idx(:,:,i,j), & ! (in)
504  interp_xih2p_coef(:, i,j), & ! (in)
505  lnpresh(:,i,j), lnpaxis(:), & ! (in)
506  var(:,i,j), & ! (in)
507  var_p(:,i,j) ) ! (out)
508  enddo
509  enddo
510 
511  return

References scale_interp::interp_interp1d().

Referenced by scale_file_history_cartesc::file_history_cartesc_truncate_3d().

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ interp_available

logical, public scale_interp_vert::interp_available = .false.

topography exists & vertical interpolation has meaning?

Definition at line 41 of file scale_interp_vert.F90.

41  logical, public :: INTERP_available = .false.

Referenced by scale_file_history_cartesc::file_history_cartesc_truncate_3d(), and interp_vert_setcoef().

scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_interp::interp_factor1d
subroutine, public interp_factor1d(KA_ref, KS_ref, KE_ref, KA, KS, KE, hgt_ref, hgt, idx_k, vfact, flag_extrap)
Definition: scale_interp.F90:252
scale_interp_vert::interp_available
logical, public interp_available
topography exists & vertical interpolation has meaning?
Definition: scale_interp_vert.F90:41
scale_interp::interp_interp1d
subroutine, public interp_interp1d(KA_ref, KS_ref, KE_ref, KA, KS, KE, idx_k, vfact, hgt_ref, hgt, val_ref, val, spline, logwgt)
Definition: scale_interp.F90:1134