SCALE-RM
scale_interp_vert.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
19  !-----------------------------------------------------------------------------
20  implicit none
21  private
22  !-----------------------------------------------------------------------------
23  !
24  !++ Public procedure
25  !
26  public :: interp_vert_setcoef
27  public :: interp_vert_xi2z
28  public :: interp_vert_xih2zh
29  public :: interp_vert_z2xi
30  public :: interp_vert_zh2xih
31 
32  public :: interp_vert_alloc_pres
33  public :: interp_vert_setcoef_pres
34  public :: interp_vert_xi2p
35  public :: interp_vert_xih2p
36 
37  !-----------------------------------------------------------------------------
38  !
39  !++ Public parameters & variables
40  !
41  logical, public :: interp_available = .false.
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private parameters & variables
50  !
51  ! full level
52  integer, private, allocatable :: interp_xi2z_idx (:,:,:,:)
53  real(rp), private, allocatable :: interp_xi2z_coef (:, :,:)
54  integer, private, allocatable :: interp_z2xi_idx (:,:,:,:)
55  real(rp), private, allocatable :: interp_z2xi_coef (:, :,:)
56 
57  integer, private, allocatable :: interp_xi2p_idx (:,:,:,:)
58  real(rp), private, allocatable :: interp_xi2p_coef (:, :,:)
59 
60  ! half level
61  integer, private, allocatable :: interp_xih2zh_idx (:,:,:,:)
62  real(rp), private, allocatable :: interp_xih2zh_coef(:, :,:)
63  integer, private, allocatable :: interp_zh2xih_idx (:,:,:,:)
64  real(rp), private, allocatable :: interp_zh2xih_coef(:, :,:)
65 
66  integer, private, allocatable :: interp_xih2p_idx (:,:,:,:)
67  real(rp), private, allocatable :: interp_xih2p_coef (:, :,:)
68 
69  ! log pressure
70  real(rp), private, allocatable :: lnpres (:,:,:)
71  real(rp), private, allocatable :: lnpresh(:,:,:)
72  real(rp), private, allocatable :: lnpaxis(:)
73 
74 
75  !-----------------------------------------------------------------------------
76 contains
77  !-----------------------------------------------------------------------------
79  subroutine interp_vert_setcoef( &
80  KA, KS, KE, &
81  IA, IS, IE, &
82  JA, JS, JE, &
83  TOPO_exist, &
84  Xi, Xih, &
85  Z, Zh )
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
175  end subroutine interp_vert_setcoef
176 
177  !-----------------------------------------------------------------------------
178  subroutine interp_vert_xi2z( &
179  KA, KS, KE, &
180  IA, IS, IE, &
181  JA, JS, JE, &
182  Xi, &
183  Z, &
184  var, &
185  var_Z )
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
216  end subroutine interp_vert_xi2z
217 
218  !-----------------------------------------------------------------------------
219  subroutine interp_vert_z2xi( &
220  KA, KS, KE, &
221  IA, IS, IE, &
222  JA, JS, JE, &
223  Z, &
224  Xi, &
225  var, &
226  var_Xi )
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
257  end subroutine interp_vert_z2xi
258 
259  !-----------------------------------------------------------------------------
260  subroutine interp_vert_xih2zh( &
261  KA, KS, KE, &
262  IA, IS, IE, &
263  JA, JS, JE, &
264  Xih, &
265  Zh, &
266  var, &
267  var_Z )
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
298  end subroutine interp_vert_xih2zh
299 
300  !-----------------------------------------------------------------------------
301  subroutine interp_vert_zh2xih( &
302  KA, KS, KE, &
303  IA, IS, IE, &
304  JA, JS, JE, &
305  Zh, &
306  Xih, &
307  var, &
308  var_Xi )
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
339  end subroutine interp_vert_zh2xih
340 
341  !-----------------------------------------------------------------------------
343  subroutine interp_vert_alloc_pres( &
344  Kpres, KA, IA, JA )
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
362  end subroutine interp_vert_alloc_pres
363 
364  !-----------------------------------------------------------------------------
365  subroutine interp_vert_setcoef_pres( &
366  Kpres, &
367  KA, KS, KE, &
368  IA, IS, IE, &
369  JA, JS, JE, &
370  PRES, &
371  PRESh, &
372  SFC_PRES, &
373  Paxis )
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
433  end subroutine interp_vert_setcoef_pres
434 
435  !-----------------------------------------------------------------------------
436  subroutine interp_vert_xi2p( &
437  Kpres, &
438  KA, KS, KE, &
439  IA, IS, IE, &
440  JA, JS, JE, &
441  var, &
442  var_P )
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
472  end subroutine interp_vert_xi2p
473 
474  !-----------------------------------------------------------------------------
475  subroutine interp_vert_xih2p( &
476  Kpres, &
477  KA, KS, KE, &
478  IA, IS, IE, &
479  JA, JS, JE, &
480  var, &
481  var_P )
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
512  end subroutine interp_vert_xih2p
513 
514 end module scale_interp_vert
scale_precision
module PRECISION
Definition: scale_precision.F90:14
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_vert_z2xi
subroutine, public interp_vert_z2xi(KA, KS, KE, IA, IS, IE, JA, JS, JE, Z, Xi, var, var_Xi)
Definition: scale_interp_vert.F90:227
scale_interp_vert::interp_vert_alloc_pres
subroutine, public interp_vert_alloc_pres(Kpres, KA, IA, JA)
Setup.
Definition: scale_interp_vert.F90:345
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_io
module STDIO
Definition: scale_io.F90:10
scale_interp_vert::interp_vert_zh2xih
subroutine, public interp_vert_zh2xih(KA, KS, KE, IA, IS, IE, JA, JS, JE, Zh, Xih, var, var_Xi)
Definition: scale_interp_vert.F90:309
scale_interp_vert::interp_available
logical, public interp_available
topography exists & vertical interpolation has meaning?
Definition: scale_interp_vert.F90:41
scale_interp_vert::interp_vert_xi2p
subroutine, public interp_vert_xi2p(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
Definition: scale_interp_vert.F90:443
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
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_interp_vert::interp_vert_xi2z
subroutine, public interp_vert_xi2z(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xi, Z, var, var_Z)
Definition: scale_interp_vert.F90:186
scale_interp_vert::interp_vert_setcoef_pres
subroutine, public interp_vert_setcoef_pres(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, PRES, PRESh, SFC_PRES, Paxis)
Definition: scale_interp_vert.F90:374
scale_interp_vert::interp_vert_xih2p
subroutine, public interp_vert_xih2p(Kpres, KA, KS, KE, IA, IS, IE, JA, JS, JE, var, var_P)
Definition: scale_interp_vert.F90:482
scale_interp_vert::interp_vert_setcoef
subroutine, public interp_vert_setcoef(KA, KS, KE, IA, IS, IE, JA, JS, JE, TOPO_exist, Xi, Xih, Z, Zh)
Setup.
Definition: scale_interp_vert.F90:86
scale_interp_vert::interp_vert_xih2zh
subroutine, public interp_vert_xih2zh(KA, KS, KE, IA, IS, IE, JA, JS, JE, Xih, Zh, var, var_Z)
Definition: scale_interp_vert.F90:268
scale_interp_vert
module INTERPOLATION vertical
Definition: scale_interp_vert.F90:11