NCEPLIBS-ip  4.1.0
ipolatev.F90
Go to the documentation of this file.
1 
4 
14  use ip_grid_mod
15 
16  implicit none
17 
18  private
20 
21  interface ipolatev
22  module procedure ipolatev_grib1
23  module procedure ipolatev_grib1_single_field
24  module procedure ipolatev_grib2
25  module procedure ipolatev_grib2_single_field
26  end interface ipolatev
27 
28 contains
29 
67  SUBROUTINE ipolatev_grid(IP,IPOPT,grid_in,grid_out, &
68  MI,MO,KM,IBI,LI,UI,VI, &
69  NO,RLAT,RLON,CROT,SROT,IBO,LO,UO,VO,IRET)
70  class(ip_grid), intent(in) :: grid_in, grid_out
71  INTEGER, INTENT(IN ) :: IP, IPOPT(20), IBI(KM)
72  INTEGER, INTENT(IN ) :: KM, MI, MO
73  INTEGER, INTENT( OUT) :: IBO(KM), IRET, NO
74  !
75  LOGICAL*1, INTENT(IN ) :: LI(MI,KM)
76  LOGICAL*1, INTENT( OUT) :: LO(MO,KM)
77  !
78  REAL, INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
79  REAL, INTENT(INOUT) :: CROT(MO),SROT(MO)
80  REAL, INTENT(INOUT) :: RLAT(MO),RLON(MO)
81  REAL, INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
82 
83  select case(ip)
84  case(bilinear_interp_id)
85  CALL interpolate_bilinear(ipopt,grid_in,grid_out, &
86  mi,mo,km,ibi,li,ui,vi,&
87  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
88  case(bicubic_interp_id)
89  CALL interpolate_bicubic(ipopt,grid_in,grid_out,mi,mo,km,ibi,li,ui,vi,&
90  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
91  case(neighbor_interp_id)
92  CALL interpolate_neighbor(ipopt,grid_in,grid_out,mi,mo,km,ibi,li,ui,vi,&
93  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
94  case(budget_interp_id)
95  CALL interpolate_budget(ipopt,grid_in,grid_out,mi,mo,km,ibi,li,ui,vi,&
96  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
97  case(spectral_interp_id)
98  CALL interpolate_spectral(ipopt,grid_in,grid_out, &
99  mi,mo,km,ibi,ui,vi,&
100  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
102  CALL interpolate_neighbor_budget(ipopt,grid_in,grid_out,mi,mo,km,ibi,li,ui,vi,&
103  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
104  case default
105  print *, "unrecognized interpolation option: ", ip
106  error stop
107  ! IF(IGDTNUMO.GE.0) NO=0
108  ! DO K=1,KM
109  ! IBO(K)=1
110  ! DO N=1,NO
111  ! LO(N,K)=.FALSE.
112  ! UO(N,K)=0.
113  ! VO(N,K)=0.
114  ! ENDDO
115  ! ENDDO
116  ! IRET=1
117  end select
118 
119  end subroutine ipolatev_grid
120 
382  subroutine ipolatev_grib2(ip,ipopt,igdtnumi,igdtmpli,igdtleni, &
383  igdtnumo,igdtmplo,igdtleno, &
384  mi,mo,km,ibi,li,ui,vi, &
385  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
386  USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool
387  INTEGER(C_INT), INTENT(IN ) :: IP, IPOPT(20), IBI(KM)
388  INTEGER(C_INT), INTENT(IN ) :: KM, MI, MO
389  INTEGER(C_INT), INTENT(IN ) :: IGDTNUMI, IGDTLENI
390  INTEGER(C_INT), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
391  INTEGER(C_INT), INTENT(IN ) :: IGDTNUMO, IGDTLENO
392  INTEGER(C_INT), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
393  INTEGER(C_INT), INTENT( OUT) :: IBO(KM), IRET, NO
394  !
395  LOGICAL(C_BOOL), INTENT(IN ) :: LI(MI,KM)
396  LOGICAL(C_BOOL), INTENT( OUT) :: LO(MO,KM)
397  !
398 #if (LSIZE==D)
399  REAL(C_DOUBLE), INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
400  REAL(C_DOUBLE), INTENT(INOUT) :: CROT(MO),SROT(MO)
401  REAL(C_DOUBLE), INTENT(INOUT) :: RLAT(MO),RLON(MO)
402  REAL(C_DOUBLE), INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
403 #elif (LSIZE==4)
404  REAL(C_FLOAT), INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
405  REAL(C_FLOAT), INTENT(INOUT) :: CROT(MO),SROT(MO)
406  REAL(C_FLOAT), INTENT(INOUT) :: RLAT(MO),RLON(MO)
407  REAL(C_FLOAT), INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
408 #endif
409  !
410 
411  type(grib2_descriptor) :: desc_in, desc_out
412  class(ip_grid), allocatable :: grid_in, grid_out
413 
414  desc_in = init_descriptor(igdtnumi, igdtleni, igdtmpli)
415  desc_out = init_descriptor(igdtnumo, igdtleno, igdtmplo)
416 
417  call init_grid(grid_in, desc_in)
418  call init_grid(grid_out, desc_out)
419 
420  CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
421  mi,mo,km,ibi,li,ui,vi,&
422  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
423 
424  end subroutine ipolatev_grib2
425 
555  subroutine ipolatev_grib1(ip,ipopt,kgdsi,kgdso,mi,mo,km,ibi,li,ui,vi, &
556  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
557  USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool
558  IMPLICIT NONE
559  !
560  INTEGER(C_INT), INTENT(IN ):: IP, IPOPT(20), IBI(KM)
561  INTEGER(C_INT), INTENT(IN ):: KM, MI, MO
562  INTEGER(C_INT), INTENT(INOUT):: KGDSI(200), KGDSO(200)
563  INTEGER(C_INT), INTENT( OUT):: IBO(KM), IRET, NO
564  !
565  LOGICAL(C_BOOL), INTENT(IN ):: LI(MI,KM)
566  LOGICAL(C_BOOL), INTENT( OUT):: LO(MO,KM)
567  !
568 #if (LSIZE==D)
569  REAL(C_DOUBLE), INTENT(IN ):: UI(MI,KM),VI(MI,KM)
570  REAL(C_DOUBLE), INTENT(INOUT):: CROT(MO),SROT(MO)
571  REAL(C_DOUBLE), INTENT(INOUT):: RLAT(MO),RLON(MO)
572  REAL(C_DOUBLE), INTENT( OUT):: UO(MO,KM),VO(MO,KM)
573 #elif (LSIZE==4)
574  REAL(C_FLOAT), INTENT(IN ):: UI(MI,KM),VI(MI,KM)
575  REAL(C_FLOAT), INTENT(INOUT):: CROT(MO),SROT(MO)
576  REAL(C_FLOAT), INTENT(INOUT):: RLAT(MO),RLON(MO)
577  REAL(C_FLOAT), INTENT( OUT):: UO(MO,KM),VO(MO,KM)
578 #endif
579  !
580  INTEGER :: KGDSI11, KGDSO11
581 
582  type(grib1_descriptor) :: desc_in, desc_out
583  class(ip_grid), allocatable :: grid_in, grid_out
584 
585  IF(kgdsi(1).EQ.203) THEN
586  kgdsi11=kgdsi(11)
587  kgdsi(11)=ior(kgdsi(11),256)
588  ENDIF
589  IF(kgdso(1).EQ.203) THEN
590  kgdso11=kgdso(11)
591  kgdso(11)=ior(kgdso(11),256)
592  ENDIF
593 
594  desc_in = init_descriptor(kgdsi)
595  desc_out = init_descriptor(kgdso)
596 
597  call init_grid(grid_in, desc_in)
598  call init_grid(grid_out, desc_out)
599 
600  CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
601  mi,mo,km,ibi,li,ui,vi,&
602  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
603 
604  IF(kgdsi(1).EQ.203) THEN
605  kgdsi(11)=kgdsi11
606  ENDIF
607  IF(kgdso(1).EQ.203) THEN
608  kgdso(11)=kgdso11
609  ENDIF
610 
611  END SUBROUTINE ipolatev_grib1
612 
663  subroutine ipolatev_grib1_single_field(ip,ipopt,kgdsi,kgdso,mi,mo,km,ibi,li,ui,vi, &
664  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
665  USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool
666  IMPLICIT NONE
667  !
668  INTEGER(C_INT), INTENT(IN ):: IP, IPOPT(20), IBI
669  INTEGER(C_INT), INTENT(IN ):: KM, MI, MO
670  INTEGER(C_INT), INTENT(INOUT):: KGDSI(200), KGDSO(200)
671  INTEGER(C_INT), INTENT( OUT):: IBO, IRET, NO
672  !
673  LOGICAL(C_BOOL), INTENT(IN ):: LI(MI)
674  LOGICAL(C_BOOL), INTENT( OUT):: LO(MO)
675  !
676 #if (LSIZE==D)
677  REAL(C_DOUBLE), INTENT(IN ):: UI(MI),VI(MI)
678  REAL(C_DOUBLE), INTENT(INOUT):: CROT(MO),SROT(MO)
679  REAL(C_DOUBLE), INTENT(INOUT):: RLAT(MO),RLON(MO)
680  REAL(C_DOUBLE), INTENT( OUT):: UO(MO),VO(MO)
681 #elif (LSIZE==4)
682  REAL(C_FLOAT), INTENT(IN ):: UI(MI),VI(MI)
683  REAL(C_FLOAT), INTENT(INOUT):: CROT(MO),SROT(MO)
684  REAL(C_FLOAT), INTENT(INOUT):: RLAT(MO),RLON(MO)
685  REAL(C_FLOAT), INTENT( OUT):: UO(MO),VO(MO)
686 #endif
687  !
688  INTEGER :: KGDSI11, KGDSO11
689 
690  type(grib1_descriptor) :: desc_in, desc_out
691  class(ip_grid), allocatable :: grid_in, grid_out
692  integer :: ibo_array(1)
693 
694  ! Can't pass expression (e.g. [ibo]) to intent(out) argument.
695  ! Initialize placeholder array of size 1 to make rank match.
696  ibo_array(1) = ibo
697 
698  IF(kgdsi(1).EQ.203) THEN
699  kgdsi11=kgdsi(11)
700  kgdsi(11)=ior(kgdsi(11),256)
701  ENDIF
702  IF(kgdso(1).EQ.203) THEN
703  kgdso11=kgdso(11)
704  kgdso(11)=ior(kgdso(11),256)
705  ENDIF
706 
707  desc_in = init_descriptor(kgdsi)
708  desc_out = init_descriptor(kgdso)
709 
710  call init_grid(grid_in, desc_in)
711  call init_grid(grid_out, desc_out)
712 
713  CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
714  mi,mo,km,[ibi],li,ui,vi,&
715  no,rlat,rlon,crot,srot,ibo_array,lo,uo,vo,iret)
716 
717  ibo = ibo_array(1)
718 
719  IF(kgdsi(1).EQ.203) THEN
720  kgdsi(11)=kgdsi11
721  ENDIF
722  IF(kgdso(1).EQ.203) THEN
723  kgdso(11)=kgdso11
724  ENDIF
725 
726  END SUBROUTINE ipolatev_grib1_single_field
727 
808  subroutine ipolatev_grib2_single_field(ip,ipopt,igdtnumi,igdtmpli,igdtleni, &
809  igdtnumo,igdtmplo,igdtleno, &
810  mi,mo,km,ibi,li,ui,vi, &
811  no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
812  USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool
813  INTEGER(C_INT), INTENT(IN ) :: IP, IPOPT(20), IBI
814  INTEGER(C_INT), INTENT(IN ) :: KM, MI, MO
815  INTEGER(C_INT), INTENT(IN ) :: IGDTNUMI, IGDTLENI
816  INTEGER(C_INT), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
817  INTEGER(C_INT), INTENT(IN ) :: IGDTNUMO, IGDTLENO
818  INTEGER(C_INT), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
819  INTEGER(C_INT), INTENT( OUT) :: IBO, IRET, NO
820  !
821  LOGICAL(C_BOOL), INTENT(IN ) :: LI(MI)
822  LOGICAL(C_BOOL), INTENT( OUT) :: LO(MO)
823  !
824 #if (LSIZE==d)
825  REAL(C_DOUBLE), INTENT(IN ) :: UI(MI),VI(MI)
826  REAL(C_DOUBLE), INTENT(INOUT) :: CROT(MO),SROT(MO)
827  REAL(C_DOUBLE), INTENT(INOUT) :: RLAT(MO),RLON(MO)
828  REAL(C_DOUBLE), INTENT( OUT) :: UO(MO),VO(MO)
829 #elif (LSIZE==4)
830  REAL(C_FLOAT), INTENT(IN ) :: UI(MI),VI(MI)
831  REAL(C_FLOAT), INTENT(INOUT) :: CROT(MO),SROT(MO)
832  REAL(C_FLOAT), INTENT(INOUT) :: RLAT(MO),RLON(MO)
833  REAL(C_FLOAT), INTENT( OUT) :: UO(MO),VO(MO)
834 #endif
835  !
836 
837  type(grib2_descriptor) :: desc_in, desc_out
838  class(ip_grid), allocatable :: grid_in, grid_out
839  integer :: ibo_array(1)
840 
841  ! Can't pass expression (e.g. [ibo]) to intent(out) argument.
842  ! Initialize placeholder array of size 1 to make rank match.
843  ibo_array(1) = ibo
844 
845  desc_in = init_descriptor(igdtnumi, igdtleni, igdtmpli)
846  desc_out = init_descriptor(igdtnumo, igdtleno, igdtmplo)
847 
848  call init_grid(grid_in, desc_in)
849  call init_grid(grid_out, desc_out)
850 
851  CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
852  mi,mo,km,[ibi],li,ui,vi,&
853  no,rlat,rlon,crot,srot,ibo_array,lo,uo,vo,iret)
854 
855  ibo = ibo_array(1)
856 
857  end subroutine ipolatev_grib2_single_field
858 
859 end module ipolatev_mod
860 
Users derived type grid descriptor objects to abstract away the raw GRIB1 and GRIB2 grid definitions.
Routines for creating an ip_grid given a Grib descriptor.
Abstract ip_grid type.
Definition: ip_grid_mod.F90:10
Top-level module to export interpolation routines and constants.
integer, parameter, public neighbor_interp_id
integer, parameter, public bilinear_interp_id
integer, parameter, public budget_interp_id
integer, parameter, public spectral_interp_id
integer, parameter, public bicubic_interp_id
integer, parameter, public neighbor_budget_interp_id
Top-level driver for vector interpolation interpolation routine ipolatev().
Definition: ipolatev.F90:10
subroutine, public ipolatev_grib2(ip, ipopt, igdtnumi, igdtmpli, igdtleni, igdtnumo, igdtmplo, igdtleno, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor.
Definition: ipolatev.F90:386
subroutine, public ipolatev_grib2_single_field(ip, ipopt, igdtnumi, igdtmpli, igdtleni, igdtnumo, igdtmplo, igdtleno, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
This subprogram interpolates vector fields from any grid to any grid given a grib2 descriptor.
Definition: ipolatev.F90:812
subroutine ipolatev_grid(IP, IPOPT, grid_in, grid_out, MI, MO, KM, IBI, LI, UI, VI, NO, RLAT, RLON, CROT, SROT, IBO, LO, UO, VO, IRET)
Interpolates vector fields between grids given ip_grid objects.
Definition: ipolatev.F90:70
subroutine, public ipolatev_grib1(ip, ipopt, kgdsi, kgdso, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
This subprogram interpolates vector field from any grid to any grid given a grib1 Grid Descriptor Sec...
Definition: ipolatev.F90:557
subroutine, public ipolatev_grib1_single_field(ip, ipopt, kgdsi, kgdso, mi, mo, km, ibi, li, ui, vi, no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret)
Special case of ipolatev_grib1 when interpolating a single field.
Definition: ipolatev.F90:665
Descriptor representing a grib1 grib descriptor section (GDS) with an integer array.
Grib-2 descriptor containing a grib2 GDT represented by an integer array.
Abstract grid that holds fields and methods common to all grids.
Definition: ip_grid_mod.F90:52