NCEPLIBS-ip 5.3.0
All Data Structures Namespaces Files Functions Variables Pages
ipolatev.F90
Go to the documentation of this file.
1!> @file
2!> @brief Top-level driver for vector interpolation routine ipolates.
3!> @author Mark Iredell, Kyle Gerheiser @date July 2021
4
5!> @brief Top-level driver for vector interpolation interpolation
6!> routine ipolatev(). The ipolatev() subprogram is overloaded with
7!> interfaces for GRIB1 and GRIB2 descriptors.
8!>
9!> @author George Gayno, Mark Iredell, Kyle Gerheiser
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
28contains
29
30 !> Interpolates vector fields between grids given ip_grid objects.
31 !>
32 !> Calls the specific interpolation routines on the generic ip_grids
33 !> created from a GRIB1 or GRIB2 descriptor.
34 !>
35 !> @param[in] ip Interpolation method.
36 !> @param[in] ipopt Interpolation options.
37 !> @param[in] grid_in Input grid.
38 !> @param[in] grid_out Output grid object created.
39 !> @param[in] mi Skip number between input grid fields if km>1 or
40 !> dimension of input grid fields if km=1.
41 !> @param[in] mo Skip number between output grid fields if km>1 or
42 !> dimension of output grid fields if km=1.
43 !> @param[in] km Number of fields to interpolate.
44 !> @param[in] ibi Input bitmap flags.
45 !> @param[in] li Input bitmaps (if respective ibi(k)=1).
46 !> @param[in] ui Input u-component fields to interpolate.
47 !> @param[in] vi Input v-component fields to interpolate.
48 !> @param[out] no Number of output points (only if kgdso(1)<0).
49 !> @param[out] rlat Output latitudes in degrees (if kgdso(1)<0).
50 !> @param[out] rlon Output longitudes in degrees (if kgdso(1)<0).
51 !> @param[inout] crot Vector rotation cosines (if igdtnumo>=0).
52 !> @param[inout] srot Vector rotation sines (if igdtnumo>=0).
53 !> @param[out] ibo Output bitmap flags.
54 !> @param[out] lo Output bitmaps (always output).
55 !> @param[out] uo Output u-component fields interpolated.
56 !> @param[out] vo Output v-component fields interpolated.
57 !> @param[out] iret Return code.
58 !> - 0 Successful interpolation.
59 !> - 1 Unrecognized interpolation method.
60 !> - 2 Unrecognized input grid or no grid overlap.
61 !> - 3 Unrecognized output grid.
62 !> - 1x Invalid bicubic method parameters.
63 !> - 3x Invalid budget method parameters.
64 !> - 4x Invalid spectral method parameters.
65 !> @date July 2021
66 !> @author Kyle Gerheiser
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)
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)
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)
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)
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)
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
121 !> This subprogram interpolates vector fields from any grid to any
122 !> grid given a grib2 descriptor.
123 !>
124 !> This is a wrapper for ipolates_grid which converts a grib1
125 !> descriptor into an ip_grid_descriptor, which is used to create
126 !> an ip_grid. Only horizontal interpolation is performed.
127 !>
128 !> The following interpolation methods are possible:
129 !> - (ip=0) bilinear
130 !> - (ip=1) bicubic
131 !> - (ip=2) neighbor
132 !> - (ip=3) budget
133 !> - (ip=4) spectral
134 !> - (ip=6) neighbor-budget
135 !>
136 !> Some of these methods have interpolation options and/or
137 !> restrictions on the input or output grids, both of which
138 !> are documented more fully in their respective subprograms.
139 !>
140 !> Input and output grids are defined by their grib 2 grid
141 !> definition template as decoded by the ncep g2 library. The
142 !> current code recognizes the following projections, where
143 !> "igdtnumi/o" is the grib 2 grid defintion template number
144 !> for the input and output grids, respectively:
145 !> - (igdtnumi/o=00) equidistant cylindrical
146 !> - (igdtnumi/o=01) rotated equidistant cylindrical. "e" and non-"e" staggered
147 !> - (igdtnumi/o=10) mercator cylindrical
148 !> - (igdtnumi/o=20) polar stereographic azimuthal
149 !> - (igdtnumi/o=30) lambert conformal conical
150 !> - (igdtnumi/o=40) gaussian cylindrical
151 !>
152 !> As an added bonus the number of output grid points
153 !> and their latitudes and longitudes are also returned.
154 !>
155 !> On the other hand, data may be interpolated to a set of station
156 !> points if "igdtnumo"<0 (or subtracted from 255 for the budget
157 !> option), in which case the number of points and
158 !> their latitudes and longitudes must be input.
159 !>
160 !> Input bitmaps will be interpolated to output bitmaps.
161 !> Output bitmaps will also be created when the output grid
162 !> extends outside of the domain of the input grid.
163 !>
164 !> The output field is set to 0 where the output bitmap is off.
165 !>
166 !> @param[in] ip Interpolation method
167 !> - ip=0 for bilinear
168 !> - ip=1 for bicubic
169 !> - ip=2 for neighbor;
170 !> - ip=3 for budget;
171 !> - ip=4 for spectral;
172 !> - ip=6 for neighbor-budget
173 !>
174 !> @param[in] ipopt Interpolation options
175 !> - ip=0: (No options)
176 !> - ip=1: Constraint option
177 !> - ip=2: (No options)
178 !> - ip=3: Number in radius, radius weights, search radius
179 !> - ip=4: Spectral shape, spectral truncation
180 !> - ip=6: Number in radius, radius weights ...)
181 !>
182 !> @param[in] igdtnumi Grid definition template number for the input grid.
183 !> Corresponds to the gfld%igdtnum component of the ncep g2 library
184 !> gridmod data structure:
185 !> - 00 - Equidistant Cylindrical
186 !> - 01 - Rotated Equidistant cylindrical. "e" and non-"e" staggered
187 !> - 10 - Mercator Cyclindrical
188 !> - 20 - Polar Stereographic azimuthal
189 !> - 30 - Lambert Conformal Conical
190 !> - 40 - Gaussian Equidistant Cyclindrical
191 !>
192 !> @param[in] igdtmpli Grid definition template array input grid.
193 !> Corresponds to the gfld%igdtmpl component of the NCEPLIBS-g2
194 !> gridmod data structure
195 !>
196 !> Section 3 Info:
197 !>
198 !> All map projections:
199 !> 1. Shape of earth, octet 15.
200 !> 2. Scale factor of spherical earth radius, octet 16.
201 !> 3. Scaled value of radius of spherical earth, octets 17-20.
202 !> 4. Scale factor of major axis of elliptical earth, octet 21.
203 !> 5. Scaled value of major axis of elliptical earth, octets 22-25.
204 !> 6. Scale factor of minor axis of elliptical earth, octet 26.
205 !> 7: Scaled value of minor axis of elliptical earth, octets 27-30.
206 !>
207 !> Equidistant Cyclindrical:
208 !> 8. Number of points along a parallel, octs 31-34.
209 !> 9. Number of points along a meridian, octs 35-38.
210 !> 10. Basic angle of initial production domain, octets 39-42.
211 !> 11. Subdivisions of basic angle, octets 43-46.
212 !> 12. Latitude of first grid point, octets 47-50.
213 !> 13. Longitude of first grid point, octets 51-54.
214 !> 14. Resolution and component flags, octet 55.
215 !> 15. Latitude of last grid point, octets 56-59.
216 !> 16. Longitude of last grid point, octets 60-63.
217 !> 17. i-direction increment, octets 64-67.
218 !> 18. j-direction increment, octets 68-71.
219 !> 19. Scanning mode, octet 72.
220 !>
221 !> Mercator Cyclindrical:
222 !> 8. Number of points along a parallel, octs 31-34.
223 !> 9. Number of points along a meridian, octs 35-38.
224 !> 10. Latitude of first point, octets 39-42.
225 !> 11. Longitude of first point, octets 43-46.
226 !> 12. Resolution and component flags, octet 47.
227 !> 13. Tangent latitude, octets 48-51.
228 !> 14. Latitude of last point, octets 52-55.
229 !> 15. Longitude of last point, octets 56-59.
230 !> 16. Scanning mode flags, octet 60.
231 !> 17. Orientation of grid, octets 61-64.
232 !> 18. Longitudinal grid length, octets 65-68.
233 !> 19. Latitudinal grid length, octets 69-72.
234 !>
235 !> Lambert Conformal Conical:
236 !> 8. Number of points along x-axis, octs 31-34.
237 !> 9. Number of points along y-axis, octs 35-38.
238 !> 10. Latitude of first point, octets 39-42.
239 !> 11. Longitude of first point, octets 43-46.
240 !> 12. Resolution of component flag, octet 47.
241 !> 13. Latitude where grid lengths specified,octets 48-51.
242 !> 14. Longitude of meridian that is parallel to y-axis, octets 52-55.
243 !> 15. x-direction grid length, octets 56-59.
244 !> 16. y-direction grid length, octets 60-63.
245 !> 17. Projection center flag, octet 64.
246 !> 18. Scanning mode, octet 65.
247 !> 19. First tangent latitude from pole, octets 66-69.
248 !> 20. Second tangent latitude from pole, octets 70-73.
249 !> 21. Latitude of south pole of projection, octets 74-77.
250 !> 22. Longitude of south pole of projection, octets 78-81.
251 !>
252 !> Gaussian Cylindrical:
253 !> 8. Number of points along a parallel, octs 31-34.
254 !> 9. Number of points along a meridian, octs 35-38.
255 !> 10. Basic angle of initial production domain, octets 39-42.
256 !> 11. Subdivisions of basic angle, octets 43-46.
257 !> 12. Latitude of first grid point, octets 47-50.
258 !> 13. Longitude of first grid point, octets 51-54.
259 !> 14. Resolution and component flags, octet 55.
260 !> 15. Latitude of last grid point, octets 56-59.
261 !> 16. Longitude of last grid point, octets 60-63.
262 !> 17. i-direction increment, octets 64-67.
263 !> 18. Number of parallels between pole and equator, octets 68-71.
264 !> 19. Scanning mode, octet 72.
265 !>
266 !> Polar Stereographic Azimuthal:
267 !> 8. Number of points along x-axis, octets 31-34.
268 !> 9. Number of points along y-axis, octets 35-38.
269 !> 10. Latitude of first grid point, octets 39-42.
270 !> 11. Longitude of first grid point, octets 43-46.
271 !> 12. Resolution and component flags, octet 47.
272 !> 13. True latitude, octets 48-51.
273 !> 14. Orientation longitude, octets 52-55.
274 !> 15. x-direction grid length, octets 56-59.
275 !> 16. y-direction grid length, octets 60-63.
276 !> 17. Projection center flag, octet 64.
277 !> 18. Scanning mode flags, octet 65.
278 !>
279 !> Rotated Equidistant Cyclindrical:
280 !> 8. Number of points along a parallel, octs 31-34.
281 !> 9. Number of points along a meridian, octs 35-38.
282 !> 10. Basic angle of initial production domain, octets 39-42.
283 !> 11. Subdivisions of basic angle, octets 43-46.
284 !> 12. Latitude of first grid point, octets 47-50.
285 !> 13. Longitude of first grid point, octets 51-54.
286 !> 14. Resolution and component flags, octet 55.
287 !> 15. Latitude of last grid point, octets 56-59.
288 !> 16. Longitude of last grid point, octets 60-63.
289 !> 17. i-direction increment, octets 64-67.
290 !> 18. j-direction increment, octets 68-71.
291 !> 19. Scanning mode, octet 72.
292 !> 20. Latitude of southern pole of projection, octets 73-76.
293 !> 21. Longitude of southern pole of projection, octets 77-80.
294 !> 22. Angle of rotation of projection, octs 81-84.
295 !>
296 !> @param[in] igdtleni Number of elements of the grid definition
297 !> template array for the input grid. Corresponds to the gfld%igdtlen
298 !> component of the ncep g2 library gridmod data structure.
299 !>
300 !> @param[in] igdtnumo Grid definition template number for the output grid.
301 !> Corresponds to the gfld%igdtnum component of the
302 !> ncep g2 library gridmod data structure.
303 !> See "igdtnumi" for specific template definitions.
304 !> Note: igdtnumo<0 means interpolate to random station points.
305 !>
306 !> @param[in] igdtmplo Grid definition template array for the output grid.
307 !> Corresponds to the gfld%igdtmpl component of the ncep g2 library
308 !> gridmod data structure.
309 !> See "igdtmpli" for definition of array elements.
310 !>
311 !> @param[in] igdtleno Number of elements of the grid definition
312 !> template array for the output grid. Corresponds to the gfld%igdtlen
313 !> component of the ncep g2 library gridmod data structure.
314 !>
315 !> @param[in] mi Skip number between input grid fields if km>1 or
316 !> dimension of input grid fields if km=1.
317 !> @param[in] mo Skip number between output grid fields if km>1
318 !> or dimension of output grid fields if km=1.
319 !> @param[in] km Number of fields to interpolate.
320 !> @param[in] ibi Input bitmap flags.
321 !> @param[in] li Input bitmaps (if respective ibi(k)=1).
322 !> @param[in] ui Input u-component fields to interpolate.
323 !> @param[in] vi Input v-component fields to interpolate.
324 !> @param[out] no Number of output points (only if kgdso(1)<0).
325 !> @param[inout] rlat Output latitudes in degrees (if kgdso(1)<0).
326 !> @param[inout] rlon Output longitudes in degrees (if kgdso(1)<0).
327 !> @param[inout] crot Vector rotation cosines (if igdtnumo>=0).
328 !> @param[inout] srot Vector rotation sines (if igdtnumo>=0).
329 !> @param[out] ibo Output bitmap flags.
330 !> @param[out] lo Output bitmaps (always output).
331 !> @param[out] uo Output u-component fields interpolated.
332 !> @param[out] vo Output v-component fields interpolated.
333 !> @param[out] iret Return code.
334 !> - 0 Successful interpolation.
335 !> - 1 Unrecognized interpolation method.
336 !> - 2 Unrecognized input grid or no grid overlap.
337 !> - 3 Unrecognized output grid.
338 !> - 1x Invalid bicubic method parameters.
339 !> - 3x Invalid budget method parameters.
340 !> - 4x Invalid spectral method parameters.
341 !>
342 !> @note Examples demonstrating relative cpu costs.
343 !> This example is interpolating 12 levels of winds
344 !> from the 360 x 181 global grid (ncep grid 3)
345 !> to the 93 x 68 hawaiian mercator grid (ncep grid 204).
346 !>
347 !> The example times are for the c90. As a reference, the cp time
348 !> for unpacking the global 12 wind fields is 0.07 seconds.
349 !>
350 !> METHOD | IP| IPOPT | CP SECONDS
351 !> --------| --|-------------| ----------
352 !> BILINEAR| 0 | | 0.05
353 !> BICUBIC | 1 | 0 | 0.16
354 !> BICUBIC | 1 | 1 | 0.17
355 !> NEIGHBOR| 2 | | 0.02
356 !> BUDGET | 3 | -1,-1 | 0.94
357 !> SPECTRAL| 4 | 0,40 | 0.31
358 !> SPECTRAL| 4 | 1,40 | 0.33
359 !> SPECTRAL| 4 | 0,-1 | 0.59
360 !> N-BUDGET| 6 | -1,-1 | 0.31
361 !>
362 !> The spectral interpolation is fast for the mercator grid.
363 !> However, for some grids the spectral interpolation is slow.
364 !>
365 !> The following example is interpolating 12 levels of winds
366 !> from the 360 x 181 global grid (ncep grid 3)
367 !> to the 93 x 65 conus lambert conformal grid (ncep grid 211).
368 !>
369 !> METHOD | IP| IPOPT |CP SECONDS
370 !> --------| --| -------------|----------
371 !> BILINEAR| 0 | | 0.05
372 !> BICUBIC | 1 | 0 | 0.15
373 !> BICUBIC | 1 | 1 | 0.16
374 !> NEIGHBOR| 2 | | 0.02
375 !> BUDGET | 3 | -1,-1 | 0.92
376 !> SPECTRAL| 4 | 0,40 | 4.51
377 !> SPECTRAL| 4 | 1,40 | 5.77
378 !> SPECTRAL| 4 | 0,-1 | 12.60
379 !> N-BUDGET| 6 | -1,-1 | 0.33
380 !>
381 !> @author Kyle Gerheiser @date July 2021
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, c_long
387#if (LSIZE==8)
388 INTEGER(C_LONG), INTENT(IN ) :: IP, IPOPT(20), IBI(KM)
389 INTEGER(C_LONG), INTENT(IN ) :: KM, MI, MO
390 INTEGER(C_LONG), INTENT(IN ) :: IGDTNUMI, IGDTLENI
391 INTEGER(C_LONG), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
392 INTEGER(C_LONG), INTENT(IN ) :: IGDTNUMO, IGDTLENO
393 INTEGER(C_LONG), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
394 INTEGER(C_LONG), INTENT( OUT) :: IBO(KM), IRET, NO
395#else
396 INTEGER(C_INT), INTENT(IN ) :: IP, IPOPT(20), IBI(KM)
397 INTEGER(C_INT), INTENT(IN ) :: KM, MI, MO
398 INTEGER(C_INT), INTENT(IN ) :: IGDTNUMI, IGDTLENI
399 INTEGER(C_INT), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
400 INTEGER(C_INT), INTENT(IN ) :: IGDTNUMO, IGDTLENO
401 INTEGER(C_INT), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
402 INTEGER(C_INT), INTENT( OUT) :: IBO(KM), IRET, NO
403#endif
404 !
405 LOGICAL(C_BOOL), INTENT(IN ) :: LI(MI,KM)
406 LOGICAL(C_BOOL), INTENT( OUT) :: LO(MO,KM)
407 !
408#if (LSIZE==4)
409 REAL(C_FLOAT), INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
410 REAL(C_FLOAT), INTENT(INOUT) :: CROT(MO),SROT(MO)
411 REAL(C_FLOAT), INTENT(INOUT) :: RLAT(MO),RLON(MO)
412 REAL(C_FLOAT), INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
413#else
414 REAL(C_DOUBLE), INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
415 REAL(C_DOUBLE), INTENT(INOUT) :: CROT(MO),SROT(MO)
416 REAL(C_DOUBLE), INTENT(INOUT) :: RLAT(MO),RLON(MO)
417 REAL(C_DOUBLE), INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
418#endif
419 !
420
421 type(grib2_descriptor) :: desc_in, desc_out
422 class(ip_grid), allocatable :: grid_in, grid_out
423
424 desc_in = init_descriptor(igdtnumi, igdtleni, igdtmpli)
425 desc_out = init_descriptor(igdtnumo, igdtleno, igdtmplo)
426
427 call init_grid(grid_in, desc_in)
428 call init_grid(grid_out, desc_out)
429
430 CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
431 mi,mo,km,ibi,li,ui,vi,&
432 no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
433
434 end subroutine ipolatev_grib2
435
436 !> @brief This subprogram interpolates vector field from any grid
437 !> to any grid given a grib1 Grid Descriptor Section.
438 !>
439 !> @details Only horizontal interpolation is performed.
440 !> The following interpolation methods are possible:
441 !> - (ip=0) bilinear
442 !> - (ip=1) bicubic
443 !> - (ip=2) neighbor
444 !> - (ip=3) budget
445 !> - (ip=4) spectral
446 !> - (ip=6) neighbor-budget
447 !>
448 !> Some of these methods have interpolation options and/or
449 !> restrictions on the input or output grids, both of which
450 !> are documented more fully in their respective subprograms.
451 !>
452 !> The grids are defined by their grid description sections
453 !> (passed in integer form as decoded by subprogram w3fi63).
454 !>
455 !> The current code recognizes the following projections:
456 !> - (kgds(1)=000) equidistant cylindrical
457 !> - (kgds(1)=001) mercator cylindrical
458 !> - (kgds(1)=003) lambert conformal conical
459 !> - (kgds(1)=004) gaussian cylindrical
460 !> - (kgds(1)=005) polar stereographic azimuthal
461 !> - (kgds(1)=203) rotated equidistant cylindrical - e-stagger
462 !> - (kgds(1)=205) rotated equidistant cylindrical - b-stagger
463 !>
464 !> Where kgds could be either input kgdsi or output kgdso.
465 !>
466 !> As an added bonus the number of output grid points
467 !> and their latitudes and longitudes are also returned.
468 !>
469 !> On the other hand, the output can be a set of station points
470 !> if kgdso(1)<0, in which case the number of points
471 !> and their latitudes and longitudes must be input.
472 !> for the budget approach, a subsection of the grid may
473 !> be output by subtracting kgdso(1) from 255 and passing
474 !> in the latitudes and longitudes of the points.
475 !> Input bitmaps will be interpolated to output bitmaps.
476 !>
477 !> Output bitmaps will also be created when the output grid
478 !> extends outside of the domain of the input grid.
479 !> the output field is set to 0 where the output bitmap is off.
480 !>
481 !> @param ip Interpolation method
482 !> - ip = BILINEAR_INTERP_ID = 0 for bilinear
483 !> - ip = BICUBIC_INTERP_ID = 1 for bicubic
484 !> - ip = NEIGHBOR_INTERP_ID = 2 for neighbor;
485 !> - ip = BUDGET_INTERP_ID = 3 for budget;
486 !> - ip = SPECTRAL_INTERP_ID = 4 for spectral;
487 !> - ip = NEIGHBOR_BUDGET_INTERP_ID = 6 for neighbor-budget
488 !> @param ipopt Interpolation options
489 !> - ip=0 (bilinear): (No options)
490 !> - ip=1 Cbicubic): constraint option
491 !> - ip=2 (neighbor): (No options)
492 !> - ip=3 (budget): Number in radius, radius weights, search radius
493 !> - ip=4 (spectral): Spectral shape, spectral truncation
494 !> - ip=6 (neighbor-budget): Number in radius, radius weights ...)
495 !> @param[in] kgdsi Input gds parameters as decoded by w3fi63.
496 !> @param[in] kgdso Output gds parameters.
497 !> @param[in] mi Skip number between input grid fields if km>1 or
498 !> dimension of input grid fields if km=1.
499 !> @param[in] mo Skip number between output grid fields if km>1 or
500 !> dimension of output grid fields if km=1.
501 !> @param[in] km Number of fields to interpolate.
502 !> @param[in] ibi Input bitmap flags.
503 !> @param[in] li Input bitmaps (if respective ibi(k)=1).
504 !> @param[in] ui Input u-component fields to interpolate.
505 !> @param[in] vi Input v-component fields to interpolate.
506 !> @param[out] no Number of output points (only if kgdso(1)<0).
507 !> @param[out] rlat Output latitudes in degrees (if kgdso(1)<0).
508 !> @param[out] rlon Output longitudes in degrees (if kgdso(1)<0).
509 !> @param[inout] crot Vector rotation cosines (if igdtnumo>=0).
510 !> @param[inout] srot Vector rotation sines (if igdtnumo>=0).
511 !> @param[out] ibo Output bitmap flags.
512 !> @param[out] lo Output bitmaps (always output).
513 !> @param[out] uo Output u-component fields interpolated.
514 !> @param[out] vo Output v-component fields interpolated.
515 !> @param[out] iret Return code.
516 !> - 0 Successful interpolation.
517 !> - 1 Unrecognized interpolation method.
518 !> - 2 Unrecognized input grid or no grid overlap.
519 !> - 3 Unrecognized output grid.
520 !> - 1x Invalid bicubic method parameters.
521 !> - 3x Invalid budget method parameters.
522 !> - 4x Invalid spectral method parameters.
523 !>
524 !> @note Examples demonstrating relative cpu costs.
525 !> This example is interpolating 12 levels of winds
526 !> from the 360 x 181 global grid (ncep grid 3)
527 !> to the 93 x 68 hawaiian mercator grid (ncep grid 204).
528 !>
529 !> The example times are for the c90. As a reference, the cp time
530 !> for unpacking the global 12 temperature fields is 0.07 seconds.
531 !>
532 !> METHOD | IP| IPOPT | CP SECONDS
533 !> --------| --|-------------| ----------
534 !> BILINEAR| 0 | | 0.05
535 !> BICUBIC | 1 | 0 | 0.16
536 !> BICUBIC | 1 | 1 | 0.17
537 !> NEIGHBOR| 2 | | 0.02
538 !> BUDGET | 3 | -1,-1 | 0.94
539 !> SPECTRAL| 4 | 0,40 | 0.31
540 !> SPECTRAL| 4 | 1,40 | 0.33
541 !> SPECTRAL| 4 | 0,-1 | 0.59
542 !> N-BUDGET| 6 | -1,-1 | 0.31
543 !>
544 !> The spectral interpolation is fast for the mercator grid.
545 !> However, for some grids the spectral interpolation is slow.
546 !>
547 !> The following example is interpolating 12 levels of temperatures
548 !> from the 360 x 181 global grid (ncep grid 3)
549 !> to the 93 x 65 conus lambert conformal grid (ncep grid 211).
550 !>
551 !> METHOD | IP| IPOPT |CP SECONDS
552 !> --------| --| -------------|----------
553 !> BILINEAR| 0 | | 0.05
554 !> BICUBIC | 1 | 0 | 0.15
555 !> BICUBIC | 1 | 1 | 0.16
556 !> NEIGHBOR| 2 | | 0.02
557 !> BUDGET | 3 | -1,-1 | 0.92
558 !> SPECTRAL| 4 | 0,40 | 4.51
559 !> SPECTRAL| 4 | 1,40 | 5.77
560 !> SPECTRAL| 4 | 0,-1 | 12.60
561 !> N-BUDGET| 6 | -1,-1 | 0.33
562 !>
563 !> @date July 2021
564 !> @author Kyle Gerheiser
565 subroutine ipolatev_grib1(ip,ipopt,kgdsi,kgdso,mi,mo,km,ibi,li,ui,vi, &
566 no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
567 USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool, c_long
568 IMPLICIT NONE
569 !
570#if (LSIZE==8)
571 INTEGER(C_LONG), INTENT(IN ):: IP, IPOPT(20), IBI(KM)
572 INTEGER(C_LONG), INTENT(IN ):: KM, MI, MO
573 INTEGER(C_LONG), INTENT(INOUT):: KGDSI(200), KGDSO(200)
574 INTEGER(C_LONG), INTENT( OUT):: IBO(KM), IRET, NO
575#else
576 INTEGER(C_INT), INTENT(IN ):: IP, IPOPT(20), IBI(KM)
577 INTEGER(C_INT), INTENT(IN ):: KM, MI, MO
578 INTEGER(C_INT), INTENT(INOUT):: KGDSI(200), KGDSO(200)
579 INTEGER(C_INT), INTENT( OUT):: IBO(KM), IRET, NO
580#endif
581 !
582 LOGICAL(C_BOOL), INTENT(IN ):: LI(MI,KM)
583 LOGICAL(C_BOOL), INTENT( OUT):: LO(MO,KM)
584 !
585#if (LSIZE==4)
586 REAL(C_FLOAT), INTENT(IN ):: UI(MI,KM),VI(MI,KM)
587 REAL(C_FLOAT), INTENT(INOUT):: CROT(MO),SROT(MO)
588 REAL(C_FLOAT), INTENT(INOUT):: RLAT(MO),RLON(MO)
589 REAL(C_FLOAT), INTENT( OUT):: UO(MO,KM),VO(MO,KM)
590#else
591 REAL(C_DOUBLE), INTENT(IN ):: UI(MI,KM),VI(MI,KM)
592 REAL(C_DOUBLE), INTENT(INOUT):: CROT(MO),SROT(MO)
593 REAL(C_DOUBLE), INTENT(INOUT):: RLAT(MO),RLON(MO)
594 REAL(C_DOUBLE), INTENT( OUT):: UO(MO,KM),VO(MO,KM)
595#endif
596 !
597 INTEGER :: KGDSI11, KGDSO11
598
599 type(grib1_descriptor) :: desc_in, desc_out
600 class(ip_grid), allocatable :: grid_in, grid_out
601
602 IF(kgdsi(1).EQ.203) THEN
603 kgdsi11=kgdsi(11)
604 kgdsi(11)=ior(kgdsi(11),256)
605 ENDIF
606 IF(kgdso(1).EQ.203) THEN
607 kgdso11=kgdso(11)
608 kgdso(11)=ior(kgdso(11),256)
609 ENDIF
610
611 desc_in = init_descriptor(kgdsi)
612 desc_out = init_descriptor(kgdso)
613
614 call init_grid(grid_in, desc_in)
615 call init_grid(grid_out, desc_out)
616
617 CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
618 mi,mo,km,ibi,li,ui,vi,&
619 no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret)
620
621 IF(kgdsi(1).EQ.203) THEN
622 kgdsi(11)=kgdsi11
623 ENDIF
624 IF(kgdso(1).EQ.203) THEN
625 kgdso(11)=kgdso11
626 ENDIF
627
628 END SUBROUTINE ipolatev_grib1
629
630 !> Special case of ipolatev_grib1 when interpolating a single field.
631 !>
632 !> Removes the km dimension of input arrays so vectors can be passed
633 !> to ibi/ibo.
634 !>
635 !> @param ip Interpolation method
636 !> - ip = BILINEAR_INTERP_ID = 0 for bilinear
637 !> - ip = BICUBIC_INTERP_ID = 1 for bicubic
638 !> - ip = NEIGHBOR_INTERP_ID = 2 for neighbor;
639 !> - ip = BUDGET_INTERP_ID = 3 for budget;
640 !> - ip = SPECTRAL_INTERP_ID = 4 for spectral;
641 !> - ip = NEIGHBOR_BUDGET_INTERP_ID = 6 for neighbor-budget
642 !> @param ipopt Interpolation options
643 !> - ip=0 (bilinear): (No options)
644 !> - ip=1 Cbicubic): constraint option
645 !> - ip=2 (neighbor): (No options)
646 !> - ip=3 (budget): Number in radius, radius weights, search radius
647 !> - ip=4 (spectral): Spectral shape, spectral truncation
648 !> - ip=6 (neighbor-budget): Number in radius, radius weights ...)
649 !> @param[in] kgdsi Input gds parameters as decoded by w3fi63.
650 !> @param[in] kgdso Output gds parameters.
651 !> @param[in] mi Skip number between input grid fields if km>1 or
652 !> dimension of input grid fields if km=1.
653 !> @param[in] mo Skip number between output grid fields if km>1 or
654 !> dimension of output grid fields if km=1.
655 !> @param[in] km Number of fields to interpolate.
656 !> @param[in] ibi Input bitmap flags.
657 !> @param[in] li Input bitmaps (if respective ibi(k)=1).
658 !> @param[in] ui Input u-component fields to interpolate.
659 !> @param[in] vi Input v-component fields to interpolate.
660 !> @param[out] no Number of output points (only if kgdso(1)<0).
661 !> @param[out] rlat Output latitudes in degrees (if kgdso(1)<0).
662 !> @param[out] rlon Output longitudes in degrees (if kgdso(1)<0).
663 !> @param[inout] crot Vector rotation cosines (if igdtnumo>=0).
664 !> @param[inout] srot Vector rotation sines (if igdtnumo>=0).
665 !> @param[out] ibo Output bitmap flags.
666 !> @param[out] lo Output bitmaps (always output).
667 !> @param[out] uo Output u-component fields interpolated.
668 !> @param[out] vo Output v-component fields interpolated.
669 !> @param[out] iret Return code.
670 !> - 0 Successful interpolation.
671 !> - 1 Unrecognized interpolation method.
672 !> - 2 Unrecognized input grid or no grid overlap.
673 !> - 3 Unrecognized output grid.
674 !> - 1x Invalid bicubic method parameters.
675 !> - 3x Invalid budget method parameters.
676 !> - 4x Invalid spectral method parameters.
677 !>
678 !> @date Jan 2022
679 !> @author Kyle Gerheiser
680 subroutine ipolatev_grib1_single_field(ip,ipopt,kgdsi,kgdso,mi,mo,km,ibi,li,ui,vi, &
681 no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
682 USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool, c_long
683 IMPLICIT NONE
684 !
685#if (LSIZE==8)
686 INTEGER(C_LONG), INTENT(IN ):: IP, IPOPT(20), IBI
687 INTEGER(C_LONG), INTENT(IN ):: KM, MI, MO
688 INTEGER(C_LONG), INTENT(INOUT):: KGDSI(200), KGDSO(200)
689 INTEGER(C_LONG), INTENT( OUT):: IBO, IRET, NO
690#else
691 INTEGER(C_INT), INTENT(IN ):: IP, IPOPT(20), IBI
692 INTEGER(C_INT), INTENT(IN ):: KM, MI, MO
693 INTEGER(C_INT), INTENT(INOUT):: KGDSI(200), KGDSO(200)
694 INTEGER(C_INT), INTENT( OUT):: IBO, IRET, NO
695#endif
696 !
697 LOGICAL(C_BOOL), INTENT(IN ):: LI(MI)
698 LOGICAL(C_BOOL), INTENT( OUT):: LO(MO)
699 !
700#if (LSIZE==4)
701 REAL(C_FLOAT), INTENT(IN ):: UI(MI),VI(MI)
702 REAL(C_FLOAT), INTENT(INOUT):: CROT(MO),SROT(MO)
703 REAL(C_FLOAT), INTENT(INOUT):: RLAT(MO),RLON(MO)
704 REAL(C_FLOAT), INTENT( OUT):: UO(MO),VO(MO)
705#else
706 REAL(C_DOUBLE), INTENT(IN ):: UI(MI),VI(MI)
707 REAL(C_DOUBLE), INTENT(INOUT):: CROT(MO),SROT(MO)
708 REAL(C_DOUBLE), INTENT(INOUT):: RLAT(MO),RLON(MO)
709 REAL(C_DOUBLE), INTENT( OUT):: UO(MO),VO(MO)
710#endif
711 !
712 INTEGER :: KGDSI11, KGDSO11
713
714 type(grib1_descriptor) :: desc_in, desc_out
715 class(ip_grid), allocatable :: grid_in, grid_out
716 integer :: ibo_array(1)
717
718 ! Can't pass expression (e.g. [ibo]) to intent(out) argument.
719 ! Initialize placeholder array of size 1 to make rank match.
720 ibo_array(1) = ibo
721
722 IF(kgdsi(1).EQ.203) THEN
723 kgdsi11=kgdsi(11)
724 kgdsi(11)=ior(kgdsi(11),256)
725 ENDIF
726 IF(kgdso(1).EQ.203) THEN
727 kgdso11=kgdso(11)
728 kgdso(11)=ior(kgdso(11),256)
729 ENDIF
730
731 desc_in = init_descriptor(kgdsi)
732 desc_out = init_descriptor(kgdso)
733
734 call init_grid(grid_in, desc_in)
735 call init_grid(grid_out, desc_out)
736
737 CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
738 mi,mo,km,[ibi],li,ui,vi,&
739 no,rlat,rlon,crot,srot,ibo_array,lo,uo,vo,iret)
740
741 ibo = ibo_array(1)
742
743 IF(kgdsi(1).EQ.203) THEN
744 kgdsi(11)=kgdsi11
745 ENDIF
746 IF(kgdso(1).EQ.203) THEN
747 kgdso(11)=kgdso11
748 ENDIF
749
750 END SUBROUTINE ipolatev_grib1_single_field
751
752 !> This subprogram interpolates vector fields from any grid to any
753 !> grid given a grib2 descriptor.
754 !>
755 !> @param[in] ip Interpolation method
756 !> - ip=0 for bilinear
757 !> - ip=1 for bicubic
758 !> - ip=2 for neighbor;
759 !> - ip=3 for budget;
760 !> - ip=4 for spectral;
761 !> - ip=6 for neighbor-budget
762 !>
763 !> @param[in] ipopt Interpolation options
764 !> - ip=0: (No options)
765 !> - ip=1: Constraint option
766 !> - ip=2: (No options)
767 !> - ip=3: Number in radius, radius weights, search radius
768 !> - ip=4: Spectral shape, spectral truncation
769 !> - ip=6: Number in radius, radius weights ...)
770 !>
771 !> @param[in] igdtnumi Grid definition template number for the input grid.
772 !> Corresponds to the gfld%igdtnum component of the ncep g2 library
773 !> gridmod data structure:
774 !> - 00 - Equidistant Cylindrical
775 !> - 01 - Rotated Equidistant cylindrical. "e" and non-"e" staggered
776 !> - 10 - Mercator Cyclindrical
777 !> - 20 - Polar Stereographic azimuthal
778 !> - 30 - Lambert Conformal Conical
779 !> - 40 - Gaussian Equidistant Cyclindrical
780 !>
781 !> @param[in] igdtmpli Grid definition template array input grid.
782 !> Corresponds to the gfld%igdtmpl component of the NCEPLIBS-g2
783 !> gridmod data structure
784 !>
785 !> @param[in] igdtleni Number of elements of the grid definition
786 !> template array for the input grid. Corresponds to the gfld%igdtlen
787 !> component of the ncep g2 library gridmod data structure.
788 !>
789 !> @param[in] igdtnumo Grid definition template number for the output grid.
790 !> Corresponds to the gfld%igdtnum component of the
791 !> ncep g2 library gridmod data structure.
792 !> See "igdtnumi" for specific template definitions.
793 !> Note: igdtnumo<0 means interpolate to random station points.
794 !>
795 !> @param[in] igdtmplo Grid definition template array for the output grid.
796 !> Corresponds to the gfld%igdtmpl component of the ncep g2 library
797 !> gridmod data structure.
798 !> See "igdtmpli" for definition of array elements.
799 !>
800 !> @param[in] igdtleno Number of elements of the grid definition
801 !> template array for the output grid. Corresponds to the gfld%igdtlen
802 !> component of the ncep g2 library gridmod data structure.
803 !>
804 !> @param[in] mi Skip number between input grid fields if km>1 or
805 !> dimension of input grid fields if km=1.
806 !> @param[in] mo Skip number between output grid fields if km>1
807 !> or dimension of output grid fields if km=1.
808 !> @param[in] km Number of fields to interpolate.
809 !> @param[in] ibi Input bitmap flags.
810 !> @param[in] li Input bitmaps (if respective ibi(k)=1).
811 !> @param[in] ui Input u-component fields to interpolate.
812 !> @param[in] vi Input v-component fields to interpolate.
813 !> @param[out] no Number of output points (only if kgdso(1)<0).
814 !> @param[inout] rlat Output latitudes in degrees (if kgdso(1)<0).
815 !> @param[inout] rlon Output longitudes in degrees (if kgdso(1)<0).
816 !> @param[inout] crot Vector rotation cosines (if igdtnumo>=0).
817 !> @param[inout] srot Vector rotation sines (if igdtnumo>=0).
818 !> @param[out] ibo Output bitmap flags.
819 !> @param[out] lo Output bitmaps (always output).
820 !> @param[out] uo Output u-component fields interpolated.
821 !> @param[out] vo Output v-component fields interpolated.
822 !> @param[out] iret Return code.
823 !> - 0 Successful interpolation.
824 !> - 1 Unrecognized interpolation method.
825 !> - 2 Unrecognized input grid or no grid overlap.
826 !> - 3 Unrecognized output grid.
827 !> - 1x Invalid bicubic method parameters.
828 !> - 3x Invalid budget method parameters.
829 !> - 4x Invalid spectral method parameters.
830 !>
831 !> @author Eric Engle @date November 2022
832 subroutine ipolatev_grib2_single_field(ip,ipopt,igdtnumi,igdtmpli,igdtleni, &
833 igdtnumo,igdtmplo,igdtleno, &
834 mi,mo,km,ibi,li,ui,vi, &
835 no,rlat,rlon,crot,srot,ibo,lo,uo,vo,iret) bind(c)
836 USE iso_c_binding, ONLY: c_int, c_float, c_double, c_bool, c_long
837#if (LSIZE==8)
838 INTEGER(C_LONG), INTENT(IN ) :: IP, IPOPT(20), IBI
839 INTEGER(C_LONG), INTENT(IN ) :: KM, MI, MO
840 INTEGER(C_LONG), INTENT(IN ) :: IGDTNUMI, IGDTLENI
841 INTEGER(C_LONG), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
842 INTEGER(C_LONG), INTENT(IN ) :: IGDTNUMO, IGDTLENO
843 INTEGER(C_LONG), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
844 INTEGER(C_LONG), INTENT( OUT) :: IBO, IRET, NO
845#else
846 INTEGER(C_INT), INTENT(IN ) :: IP, IPOPT(20), IBI
847 INTEGER(C_INT), INTENT(IN ) :: KM, MI, MO
848 INTEGER(C_INT), INTENT(IN ) :: IGDTNUMI, IGDTLENI
849 INTEGER(C_INT), INTENT(IN ) :: IGDTMPLI(IGDTLENI)
850 INTEGER(C_INT), INTENT(IN ) :: IGDTNUMO, IGDTLENO
851 INTEGER(C_INT), INTENT(IN ) :: IGDTMPLO(IGDTLENO)
852 INTEGER(C_INT), INTENT( OUT) :: IBO, IRET, NO
853#endif
854 !
855 LOGICAL(C_BOOL), INTENT(IN ) :: LI(MI)
856 LOGICAL(C_BOOL), INTENT( OUT) :: LO(MO)
857 !
858#if (LSIZE==4)
859 REAL(C_FLOAT), INTENT(IN ) :: UI(MI),VI(MI)
860 REAL(C_FLOAT), INTENT(INOUT) :: CROT(MO),SROT(MO)
861 REAL(C_FLOAT), INTENT(INOUT) :: RLAT(MO),RLON(MO)
862 REAL(C_FLOAT), INTENT( OUT) :: UO(MO),VO(MO)
863#else
864 REAL(C_DOUBLE), INTENT(IN ) :: UI(MI),VI(MI)
865 REAL(C_DOUBLE), INTENT(INOUT) :: CROT(MO),SROT(MO)
866 REAL(C_DOUBLE), INTENT(INOUT) :: RLAT(MO),RLON(MO)
867 REAL(C_DOUBLE), INTENT( OUT) :: UO(MO),VO(MO)
868#endif
869 !
870
871 type(grib2_descriptor) :: desc_in, desc_out
872 class(ip_grid), allocatable :: grid_in, grid_out
873 integer :: ibo_array(1)
874
875 ! Can't pass expression (e.g. [ibo]) to intent(out) argument.
876 ! Initialize placeholder array of size 1 to make rank match.
877 ibo_array(1) = ibo
878
879 desc_in = init_descriptor(igdtnumi, igdtleni, igdtmpli)
880 desc_out = init_descriptor(igdtnumo, igdtleno, igdtmplo)
881
882 call init_grid(grid_in, desc_in)
883 call init_grid(grid_out, desc_out)
884
885 CALL ipolatev_grid(ip,ipopt,grid_in,grid_out, &
886 mi,mo,km,[ibi],li,ui,vi,&
887 no,rlat,rlon,crot,srot,ibo_array,lo,uo,vo,iret)
888
889 ibo = ibo_array(1)
890
891 end subroutine ipolatev_grib2_single_field
892
893end module ipolatev_mod
894
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.
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 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_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:836
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:567
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:682
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.