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_dealloc_pres
34  public :: interp_vert_setcoef_pres
35  public :: interp_vert_xi2p
36  public :: interp_vert_xih2p
37  public :: interp_vert_finalize
38 
39  !-----------------------------------------------------------------------------
40  !
41  !++ Public parameters & variables
42  !
43  logical, public :: interp_available = .false.
44 
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  !-----------------------------------------------------------------------------
50  !
51  !++ Private parameters & variables
52  !
53  ! full level
54  integer, private, allocatable :: interp_xi2z_idx (:,:,:,:)
55  real(rp), private, allocatable :: interp_xi2z_coef (:, :,:)
56  integer, private, allocatable :: interp_z2xi_idx (:,:,:,:)
57  real(rp), private, allocatable :: interp_z2xi_coef (:, :,:)
58 
59  integer, private, allocatable :: interp_xi2p_idx (:,:,:,:)
60  real(rp), private, allocatable :: interp_xi2p_coef (:, :,:)
61 
62  ! half level
63  integer, private, allocatable :: interp_xih2zh_idx (:,:,:,:)
64  real(rp), private, allocatable :: interp_xih2zh_coef(:, :,:)
65  integer, private, allocatable :: interp_zh2xih_idx (:,:,:,:)
66  real(rp), private, allocatable :: interp_zh2xih_coef(:, :,:)
67 
68  integer, private, allocatable :: interp_xih2p_idx (:,:,:,:)
69  real(rp), private, allocatable :: interp_xih2p_coef (:, :,:)
70 
71  ! log pressure
72  real(rp), private, allocatable :: lnpres (:,:,:)
73  real(rp), private, allocatable :: lnpresh(:,:,:)
74  real(rp), private, allocatable :: lnpaxis(:)
75 
76 
77  !-----------------------------------------------------------------------------
78 contains
79  !-----------------------------------------------------------------------------
81  subroutine interp_vert_setcoef( &
82  KA, KS, KE, &
83  IA, IS, IE, &
84  JA, JS, JE, &
85  TOPO_exist, &
86  Xi, Xih, &
87  Z, Zh )
88  use scale_interp, only: &
90  implicit none
91  integer, intent(in) :: ka, ks, ke
92  integer, intent(in) :: ia, is, ie
93  integer, intent(in) :: ja, js, je
94  logical, intent(in) :: topo_exist
95  real(rp), intent(in) :: xi ( ka)
96  real(rp), intent(in) :: xih(0:ka)
97  real(rp), intent(in) :: z ( ka,ia,ja)
98  real(rp), intent(in) :: zh (0:ka,ia,ja)
99 
100  integer :: i, j
101  !---------------------------------------------------------------------------
102 
103  log_newline
104  log_info("INTERP_VERT_setcoef",*) 'Setup'
105  log_info("INTERP_VERT_setcoef",*) 'No namelists.'
106 
107  interp_available = topo_exist
108 
109  log_newline
110  log_info("INTERP_VERT_setcoef",*) 'Topography exists & interpolation has meaning? : ', interp_available
111 
112  ! full level
113 
114  allocate( interp_xi2z_idx(ka,2,ia,ja) )
115  allocate( interp_xi2z_coef(ka, ia,ja) )
116  allocate( interp_z2xi_idx(ka,2,ia,ja) )
117  allocate( interp_z2xi_coef(ka, ia,ja) )
118  !$acc enter data create(INTERP_xi2z_idx, INTERP_xi2z_coef, INTERP_z2xi_idx, INTERP_z2xi_coef)
119 
120  !$acc data copyin(Z, Xi)
121 
122  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
123  !$omp shared(JA,IA,KA,KS,KE,Xi,Z,INTERP_xi2z_idx,INTERP_xi2z_coef)
124  !$acc kernels
125  !$acc loop independent
126  do j = 1, ja
127  !$acc loop independent
128  do i = 1, ia
129  call interp_factor1d( ka, ks, ke, ka, ks, ke, &
130  z(:,i,j), xi(:), & ! (in)
131  interp_xi2z_idx(:,:,i,j), & ! (out)
132  interp_xi2z_coef(:, i,j), & ! (out)
133  flag_extrap = .true. ) ! (in)
134  enddo
135  enddo
136  !$acc end kernels
137 
138  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
139  !$omp shared(JA,IA,KA,KS,KE,Z,Xi,INTERP_z2xi_idx,INTERP_z2xi_coef)
140  !$acc kernels
141  !$acc loop independent
142  do j = 1, ja
143  !$acc loop independent
144  do i = 1, ia
145  call interp_factor1d( ka, ks, ke, ka, ks, ke, &
146  xi(:), z(:,i,j), & ! (in)
147  interp_z2xi_idx(:,:,i,j), & ! (out)
148  interp_z2xi_coef(:, i,j), & ! (out)
149  flag_extrap = .true. ) ! (in)
150  enddo
151  enddo
152  !$acc end kernels
153 
154  !$acc end data
155 
156  ! half level
157 
158  allocate( interp_xih2zh_idx(ka,2,ia,ja) )
159  allocate( interp_xih2zh_coef(ka, ia,ja) )
160  allocate( interp_zh2xih_idx(ka,2,ia,ja) )
161  allocate( interp_zh2xih_coef(ka, ia,ja) )
162  !$acc enter data create(INTERP_xih2zh_idx, INTERP_xih2zh_coef, INTERP_zh2xih_idx, INTERP_zh2xih_coef)
163 
164  !$acc data copyin(Zh, Xih)
165 
166  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
167  !$omp shared(JA,IA,KA,KS,KE,Xih,Zh,INTERP_xih2zh_idx,INTERP_xih2zh_coef)
168  !$acc kernels
169  !$acc loop independent
170  do j = 1, ja
171  !$acc loop independent
172  do i = 1, ia
173  call interp_factor1d( ka, ks-1, ke, ka, ks-1, ke, &
174  zh(1:,i,j), xih(1:), & ! (in)
175  interp_xih2zh_idx(:,:,i,j), & ! (out)
176  interp_xih2zh_coef(:, i,j), & ! (out)
177  flag_extrap = .false. ) ! (in)
178  enddo
179  enddo
180  !$acc end kernels
181 
182  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
183  !$omp shared(JA,IA,KA,KS,KE,Xih,Zh,INTERP_zh2xih_idx,INTERP_zh2xih_coef)
184  !$acc kernels
185  !$acc loop independent
186  do j = 1, ja
187  !$acc loop independent
188  do i = 1, ia
189  call interp_factor1d( ka, ks-1, ke, ka, ks-1, ke, &
190  xih(1:), zh(1:,i,j), & ! (in)
191  interp_zh2xih_idx(:,:,i,j), & ! (out)
192  interp_zh2xih_coef(:, i,j), & ! (out)
193  flag_extrap = .true. ) ! (in)
194  enddo
195  enddo
196  !$acc end kernels
197 
198  !$acc end data
199 
200  return
201  end subroutine interp_vert_setcoef
202 
203  !-----------------------------------------------------------------------------
204  subroutine interp_vert_xi2z( &
205  KA, KS, KE, &
206  IA, IS, IE, &
207  JA, JS, JE, &
208  Xi, &
209  Z, &
210  var, &
211  var_Z )
212  use scale_interp, only: &
214  implicit none
215  integer, intent(in) :: ka, ks, ke
216  integer, intent(in) :: ia, is, ie
217  integer, intent(in) :: ja, js, je
218  real(rp), intent(in) :: xi (ka)
219  real(rp), intent(in) :: z (ka,ia,ja)
220  real(rp), intent(in) :: var (ka,ia,ja)
221  real(rp), intent(out) :: var_z(ka,ia,ja)
222 
223 #ifdef _OPENACC
224  real(rp) :: workr(ka,7)
225  integer :: worki(ka,2)
226 #endif
227 
228  integer :: i, j
229  !---------------------------------------------------------------------------
230 
231  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
232  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
233  !$omp shared(Xi,Z,var,var_Z) &
234  !$omp shared(INTERP_xi2z_idx,INTERP_xi2z_coef)
235  !$acc kernels copyin(Xi, Z, var) copyout(var_Z)
236  !$acc loop independent
237  do j = js, je
238  !$acc loop independent private(workr,worki)
239  do i = is, ie
240  call interp_interp1d( ka, ks, ke, ka, ks, ke, &
241 #ifdef _OPENACC
242  workr(:,:), worki(:,:), &
243 #endif
244  interp_xi2z_idx(:,:,i,j), & ! (in)
245  interp_xi2z_coef(:, i,j), & ! (in)
246  z(:,i,j), xi(:), & ! (in)
247  var(:,i,j), & ! (in)
248  var_z(:,i,j) ) ! (out)
249  end do
250  end do
251  !$acc end kernels
252 
253  return
254  end subroutine interp_vert_xi2z
255 
256  !-----------------------------------------------------------------------------
257  subroutine interp_vert_z2xi( &
258  KA, KS, KE, &
259  IA, IS, IE, &
260  JA, JS, JE, &
261  Z, &
262  Xi, &
263  var, &
264  var_Xi )
265  use scale_interp, only: &
267  implicit none
268  integer, intent(in) :: ka, ks, ke
269  integer, intent(in) :: ia, is, ie
270  integer, intent(in) :: ja, js, je
271  real(rp), intent(in) :: z (ka,ia,ja)
272  real(rp), intent(in) :: xi (ka)
273  real(rp), intent(in) :: var (ka,ia,ja)
274  real(rp), intent(out) :: var_xi(ka,ia,ja)
275 
276 #ifdef _OPENACC
277  real(rp) :: workr(ka,7)
278  integer :: worki(ka,2)
279 #endif
280 
281  integer :: i, j
282  !---------------------------------------------------------------------------
283 
284  call prof_rapstart('INTERP_interp',3)
285 
286  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
287  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
288  !$omp shared(Z,Xi,var,var_Xi) &
289  !$omp shared(INTERP_z2xi_idx,INTERP_z2xi_coef)
290  !$acc kernels copyin(Z, Xi, var) copyout(var_Xi)
291  !$acc loop independent
292  do j = js, je
293  !$acc loop independent private(workr,worki)
294  do i = is, ie
295  call interp_interp1d( ka, ks, ke, ka, ks, ke, &
296 #ifdef _OPENACC
297  workr(:,:), worki(:,:), &
298 #endif
299  interp_z2xi_idx(:,:,i,j), & ! (in)
300  interp_z2xi_coef(:, i,j), & ! (in)
301  xi(:), z(:,i,j), & ! (in)
302  var(:,i,j), & ! (in)
303  var_xi(:,i,j) ) ! (out)
304  end do
305  end do
306  !$acc end kernels
307 
308  call prof_rapend ('INTERP_interp',3)
309 
310  return
311  end subroutine interp_vert_z2xi
312 
313  !-----------------------------------------------------------------------------
314  subroutine interp_vert_xih2zh( &
315  KA, KS, KE, &
316  IA, IS, IE, &
317  JA, JS, JE, &
318  Xih, &
319  Zh, &
320  var, &
321  var_Z )
322  use scale_interp, only: &
324  implicit none
325  integer, intent(in) :: ka, ks, ke
326  integer, intent(in) :: ia, is, ie
327  integer, intent(in) :: ja, js, je
328  real(rp), intent(in) :: xih (0:ka)
329  real(rp), intent(in) :: zh (0:ka,ia,ja)
330  real(rp), intent(in) :: var (ka,ia,ja)
331  real(rp), intent(out) :: var_z(ka,ia,ja)
332 
333  real(rp) :: zht(ka)
334  real(rp) :: xiht(ka)
335 #ifdef _OPENACC
336  real(rp) :: workr(ka,7)
337  integer :: worki(ka,2)
338 #endif
339 
340  integer :: i, j, k
341  !---------------------------------------------------------------------------
342 
343  !$acc data copyin(Zh, var) copyout(var_Z) create(Xiht)
344 
345  !$kernels
346  do k = ks-1, ke
347  xiht(k) = xih(k)
348  end do
349  !$end kernels
350 
351  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
352  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
353  !$omp shared(Xih,Xiht,Zh,var,var_Z) &
354  !$omp shared(INTERP_xih2zh_idx,INTERP_xih2zh_coef) &
355  !$omp private(Zht)
356  !$acc kernels
357  !$acc loop independent
358  do j = js, je
359  !$acc loop independent private(workr,worki,Zht)
360  do i = is, ie
361  do k = ks-1, ke
362  zht(k) = zh(k,i,j)
363  end do
364  call interp_interp1d( ka, ks-1, ke, ka, ks-1, ke, &
365 #ifdef _OPENACC
366  workr(:,:), worki(:,:), &
367 #endif
368  interp_xih2zh_idx(:,:,i,j), & ! (in)
369  interp_xih2zh_coef(:, i,j), & ! (in)
370  zht(:), xiht(:), & ! (in)
371  var(:,i,j), & ! (in)
372  var_z(:,i,j) ) ! (out)
373  enddo
374  enddo
375  !$acc end kernels
376 
377  !$acc end data
378 
379  return
380  end subroutine interp_vert_xih2zh
381 
382  !-----------------------------------------------------------------------------
383  subroutine interp_vert_zh2xih( &
384  KA, KS, KE, &
385  IA, IS, IE, &
386  JA, JS, JE, &
387  Zh, &
388  Xih, &
389  var, &
390  var_Xi )
391  use scale_interp, only: &
393  implicit none
394  integer, intent(in) :: ka, ks, ke
395  integer, intent(in) :: ia, is, ie
396  integer, intent(in) :: ja, js, je
397  real(rp), intent(in) :: zh (0:ka,ia,ja)
398  real(rp), intent(in) :: xih (0:ka)
399  real(rp), intent(in) :: var (ka,ia,ja)
400  real(rp), intent(out) :: var_xi(ka,ia,ja)
401 
402  real(rp) :: zht(ka)
403  real(rp) :: xiht(ka)
404 #ifdef _OPENACC
405  real(rp) :: workr(ka,7)
406  integer :: worki(ka,2)
407 #endif
408 
409  integer :: k, i, j
410  !---------------------------------------------------------------------------
411 
412  !$acc data copyin(Zh, Xih, var) copyout(var_Xi) create(Xih)
413 
414  !$kernels
415  do k = ks-1, ke
416  xiht(k) = xih(k)
417  end do
418  !$end kernels
419 
420  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
421  !$omp shared(KA,KS,KE,IS,IE,JS,JE) &
422  !$omp shared(Zh,Xih,Xiht,var,var_Xi) &
423  !$omp shared(INTERP_zh2xih_idx,INTERP_zh2xih_coef) &
424  !$omp private(Zht)
425  !$acc kernels
426  !$acc loop independent
427  do j = js, je
428  !$acc loop independent private(workr,worki,Zht)
429  do i = is, ie
430  do k = ks-1, ke
431  zht(k) = zh(k,i,j)
432  end do
433  call interp_interp1d( ka, ks-1, ke, ka, ks-1, ke, &
434 #ifdef _OPENACC
435  workr(:,:), worki(:,:), &
436 #endif
437  interp_zh2xih_idx(:,:,i,j), & ! (in)
438  interp_zh2xih_coef(:, i,j), & ! (in)
439  xiht(:), zht(:), & ! (in)
440  var(:,i,j), & ! (in)
441  var_xi(:,i,j) ) ! (out)
442  enddo
443  enddo
444  !$acc end kernels
445 
446  !$acc end data
447 
448  return
449  end subroutine interp_vert_zh2xih
450 
451  !-----------------------------------------------------------------------------
453  subroutine interp_vert_alloc_pres( &
454  Kpres, KA, IA, JA )
455  implicit none
456  integer, intent(in) :: kpres
457  integer, intent(in) :: ka
458  integer, intent(in) :: ia
459  integer, intent(in) :: ja
460  !---------------------------------------------------------------------------
461 
462  allocate( interp_xi2p_idx(kpres,2,ia,ja) )
463  allocate( interp_xi2p_coef(kpres, ia,ja) )
464  allocate( interp_xih2p_idx(kpres,2,ia,ja) )
465  allocate( interp_xih2p_coef(kpres, ia,ja) )
466  !$acc enter data create(INTERP_xi2p_idx, INTERP_xi2p_coef, INTERP_xih2p_idx, INTERP_xih2p_coef)
467 
468  allocate( lnpres(ka,ia,ja) )
469  allocate( lnpresh(ka,ia,ja) )
470  allocate( lnpaxis(kpres) )
471  !$acc enter data create(LnPRES, LnPRESh, LnPaxis)
472 
473  return
474  end subroutine interp_vert_alloc_pres
475 
476  !-----------------------------------------------------------------------------
478  subroutine interp_vert_dealloc_pres
479  implicit none
480  !---------------------------------------------------------------------------
481 
482  !$acc exit data delete(INTERP_xi2p_idx, INTERP_xi2p_coef, INTERP_xih2p_idx, INTERP_xih2p_coef)
483  deallocate( interp_xi2p_idx )
484  deallocate( interp_xi2p_coef )
485  deallocate( interp_xih2p_idx )
486  deallocate( interp_xih2p_coef )
487 
488  !$acc exit data delete(LnPRES, LnPRESh, LnPaxis)
489  deallocate( lnpres )
490  deallocate( lnpresh )
491  deallocate( lnpaxis )
492 
493  return
494  end subroutine interp_vert_dealloc_pres
495 
496  !-----------------------------------------------------------------------------
497  subroutine interp_vert_setcoef_pres( &
498  Kpres, &
499  KA, KS, KE, &
500  IA, IS, IE, &
501  JA, JS, JE, &
502  PRES, &
503  PRESh, &
504  SFC_PRES, &
505  Paxis )
506  use scale_interp, only: &
508  implicit none
509 
510  integer, intent(in) :: kpres
511  integer, intent(in) :: ka, ks, ke
512  integer, intent(in) :: ia, is, ie
513  integer, intent(in) :: ja, js, je
514  real(rp), intent(in) :: pres ( ka,ia,ja) ! pressure in Xi coordinate [Pa]
515  real(rp), intent(in) :: presh (0:ka,ia,ja) ! pressure in Xi coordinate [Pa], layer interface
516  real(rp), intent(in) :: sfc_pres( ia,ja) ! surface pressure [Pa]
517  real(rp), intent(in) :: paxis (kpres) ! pressure level to output [Pa]
518 
519  integer :: k, i, j
520  !---------------------------------------------------------------------------
521 
522  !$acc data copyin(PRES, PRESh, SFC_PRES, Paxis)
523 
524  ! full level
525 
526 !OCL SERIAL
527  !$acc kernels
528  do k = 1, kpres
529  lnpaxis(k) = - log( paxis(k) )
530  end do
531  !$acc end kernels
532 
533  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
534  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
535  !$omp shared(PRES,LnPRES,LnPaxis) &
536  !$omp shared(INTERP_xi2p_idx,INTERP_xi2p_coef)
537  !$acc kernels
538  !$acc loop independent
539  do j = js, je
540  !$acc loop independent
541  do i = is, ie
542  do k = ks, ke
543  lnpres(k,i,j) = - log( pres(k,i,j) )
544  end do
545  call interp_factor1d( ka, ks, ke, kpres, 1, kpres, &
546  lnpres(:,i,j), lnpaxis(:), & ! (in)
547  interp_xi2p_idx(:,:,i,j), & ! (out)
548  interp_xi2p_coef(:, i,j), & ! (out)
549  flag_extrap = .false. ) ! (in)
550  enddo
551  enddo
552  !$acc end kernels
553 
554  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
555  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
556  !$omp shared(PRESh,LnPRESh,LnPaxis) &
557  !$omp shared(INTERP_xih2p_idx,INTERP_xih2p_coef)
558  !$acc kernels
559  !$acc loop independent
560  do j = js, je
561  !$acc loop independent
562  do i = is, ie
563  do k = ks-1, ke
564  lnpresh(k,i,j) = - log( presh(k,i,j) )
565  end do
566  call interp_factor1d( ka, ks-1, ke, kpres, 1, kpres, &
567  lnpresh(:,i,j), lnpaxis(:), & ! (in)
568  interp_xih2p_idx(:,:,i,j), & ! (out)
569  interp_xih2p_coef(:, i,j), & ! (out)
570  flag_extrap = .false. ) ! (in)
571  enddo
572  enddo
573  !$acc end kernels
574 
575  !$acc end data
576 
577  return
578  end subroutine interp_vert_setcoef_pres
579 
580  !-----------------------------------------------------------------------------
581  subroutine interp_vert_xi2p( &
582  Kpres, &
583  KA, KS, KE, &
584  IA, IS, IE, &
585  JA, JS, JE, &
586  var, &
587  var_P )
588  use scale_interp, only: &
590  implicit none
591  integer, intent(in) :: kpres
592  integer, intent(in) :: ka, ks, ke
593  integer, intent(in) :: ia, is, ie
594  integer, intent(in) :: ja, js, je
595  real(rp), intent(in) :: var (ka ,ia,ja)
596  real(rp), intent(out) :: var_p(kpres,ia,ja)
597 
598 #ifdef _OPENACC
599  real(rp) :: workr(ka,7)
600  integer :: worki(ka,2)
601 #endif
602 
603  integer :: i, j
604  !---------------------------------------------------------------------------
605 
606  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
607  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
608  !$omp shared(LnPRES,LnPaxis,var,var_P) &
609  !$omp shared(INTERP_xi2p_idx,INTERP_xi2p_coef)
610  !$acc kernels copyin(var) copyout(var_P)
611  !$acc loop independent
612  do j = js, je
613  !$acc loop independent private(workr,worki)
614  do i = is, ie
615  call interp_interp1d( ka, ks, ke, kpres, 1, kpres, &
616 #ifdef _OPENACC
617  workr(:,:), worki(:,:), &
618 #endif
619  interp_xi2p_idx(:,:,i,j), & ! (in)
620  interp_xi2p_coef(:, i,j), & ! (in)
621  lnpres(:,i,j), lnpaxis(:), & ! (in)
622  var(:,i,j), & ! (in)
623  var_p(:,i,j) ) ! (out)
624  enddo
625  enddo
626  !$acc end kernels
627 
628  return
629  end subroutine interp_vert_xi2p
630 
631  !-----------------------------------------------------------------------------
632  subroutine interp_vert_xih2p( &
633  Kpres, &
634  KA, KS, KE, &
635  IA, IS, IE, &
636  JA, JS, JE, &
637  var, &
638  var_P )
639  use scale_interp, only: &
641  implicit none
642 
643  integer, intent(in) :: kpres
644  integer, intent(in) :: ka, ks, ke
645  integer, intent(in) :: ia, is, ie
646  integer, intent(in) :: ja, js, je
647  real(rp), intent(in) :: var (ka ,ia,ja)
648  real(rp), intent(out) :: var_p(kpres,ia,ja)
649 
650 #ifdef _OPENACC
651  real(rp) :: workr(ka,7)
652  integer :: worki(ka,2)
653 #endif
654 
655  integer :: k, i, j
656  !---------------------------------------------------------------------------
657 
658  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
659  !$omp shared(KA,KS,KE,Kpres,IS,IE,JS,JE) &
660  !$omp shared(LnPRESh,LnPaxis,var,var_P) &
661  !$omp shared(INTERP_xih2p_idx,INTERP_xih2p_coef)
662  !$acc kernels copyin(var) copyout(var_P)
663  !$acc loop independent
664  do j = js, je
665  !$acc loop independent private(workr,worki)
666  do i = is, ie
667  call interp_interp1d( ka, ks-1, ke, kpres, 1, kpres, &
668 #ifdef _OPENACC
669  workr(:,:), worki(:,:), &
670 #endif
671  interp_xih2p_idx(:,:,i,j), & ! (in)
672  interp_xih2p_coef(:, i,j), & ! (in)
673  lnpresh(:,i,j), lnpaxis(:), & ! (in)
674  var(:,i,j), & ! (in)
675  var_p(:,i,j) ) ! (out)
676  enddo
677  enddo
678  !$acc end kernels
679 
680  return
681  end subroutine interp_vert_xih2p
682 
683  !-----------------------------------------------------------------------------
685  subroutine interp_vert_finalize
686  implicit none
687  !---------------------------------------------------------------------------
688 
689  !$acc exit data delete(INTERP_xi2z_idx, INTERP_xi2z_coef, INTERP_z2xi_idx, INTERP_z2xi_coef)
690  deallocate( interp_xi2z_idx )
691  deallocate( interp_xi2z_coef )
692  deallocate( interp_z2xi_idx )
693  deallocate( interp_z2xi_coef )
694 
695  !$acc exit data delete(INTERP_xih2zh_idx, INTERP_xih2zh_coef, INTERP_zh2xih_idx, INTERP_zh2xih_coef)
696  deallocate( interp_xih2zh_idx )
697  deallocate( interp_xih2zh_coef )
698  deallocate( interp_zh2xih_idx )
699  deallocate( interp_zh2xih_coef )
700 
701  return
702  end subroutine interp_vert_finalize
703 
704 end module scale_interp_vert
scale_interp_vert::interp_vert_finalize
subroutine, public interp_vert_finalize
Finalize.
Definition: scale_interp_vert.F90:686
scale_interp_vert::interp_vert_dealloc_pres
subroutine, public interp_vert_dealloc_pres
Finalize for pressure coordinate.
Definition: scale_interp_vert.F90:479
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_prof::prof_rapstart
subroutine, public prof_rapstart(rapname_base, level, disable_barrier)
Start raptime.
Definition: scale_prof.F90:174
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:256
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:265
scale_interp_vert::interp_vert_alloc_pres
subroutine, public interp_vert_alloc_pres(Kpres, KA, IA, JA)
Setup for pressure coordinate.
Definition: scale_interp_vert.F90:455
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:391
scale_interp_vert::interp_available
logical, public interp_available
topography exists & vertical interpolation has meaning?
Definition: scale_interp_vert.F90:43
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:588
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:1263
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:212
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:506
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:639
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:88
scale_prof::prof_rapend
subroutine, public prof_rapend(rapname_base, level, disable_barrier)
Save raptime.
Definition: scale_prof.F90:246
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:322
scale_interp_vert
module INTERPOLATION vertical
Definition: scale_interp_vert.F90:11