35 REAL :: TINYREAL=tiny(1.0)
101 MI,MO,KM,IBI,LI,GI, &
102 NO,RLAT,RLON,IBO,LO,GO,IRET)
103 class(ip_grid),
intent(in) :: grid_in, grid_out
104 INTEGER,
INTENT(IN ) :: IPOPT(20)
105 INTEGER,
INTENT(IN ) :: MI,MO,KM
106 INTEGER,
INTENT(IN ) :: IBI(KM)
107 INTEGER,
INTENT(INOUT) :: NO
108 INTEGER,
INTENT( OUT) :: IRET, IBO(KM)
110 LOGICAL*1,
INTENT(IN ) :: LI(MI,KM)
111 LOGICAL*1,
INTENT( OUT) :: LO(MO,KM)
113 REAL,
INTENT(IN ) :: GI(MI,KM)
114 REAL,
INTENT(INOUT) :: RLAT(MO),RLON(MO)
115 REAL,
INTENT( OUT) :: GO(MO,KM)
117 REAL,
PARAMETER :: FILL=-9999.
119 INTEGER :: I1,J1,IXS,JXS
120 INTEGER :: MSPIRAL,N,K,NK
122 INTEGER :: MX,KXS,KXT,IX,JX,NX
124 LOGICAL :: SAME_GRIDI, SAME_GRIDO
126 REAL :: XPTS(MO),YPTS(MO)
127 logical :: to_station_points
129 INTEGER,
SAVE :: NOX=-1,iretx=-1
130 INTEGER,
ALLOCATABLE,
SAVE :: NXY(:)
131 REAL,
ALLOCATABLE,
SAVE :: RLATX(:),RLONX(:),XPTSX(:),YPTSX(:)
132 class(ip_grid),
allocatable,
save :: prev_grid_in, prev_grid_out
136 mspiral=max(ipopt(1),1)
138 if (.not.
allocated(prev_grid_in) .or. .not.
allocated(prev_grid_out))
then
139 allocate(prev_grid_in, source = grid_in)
140 allocate(prev_grid_out, source = grid_out)
145 same_gridi = grid_in == prev_grid_in
146 same_grido = grid_out == prev_grid_out
148 if (.not. same_gridi .or. .not. same_grido)
then
149 deallocate(prev_grid_in)
150 deallocate(prev_grid_out)
152 allocate(prev_grid_in, source = grid_in)
153 allocate(prev_grid_out, source = grid_out)
157 select type(grid_out)
158 type is(ip_station_points_grid)
159 to_station_points = .true.
161 to_station_points = .false.
165 IF(iret.EQ.0.AND.(to_station_points.OR..NOT.same_gridi.OR..NOT.same_grido))
THEN
168 CALL gdswzd(grid_out, 0,mo,fill,xpts,ypts,rlon,rlat,no)
172 CALL gdswzd(grid_in,-1,no,fill,xpts,ypts,rlon,rlat,nv)
173 IF(iret.EQ.0.AND.nv.EQ.0) iret=2
177 IF(nox.GE.0)
DEALLOCATE(rlatx,rlonx,xptsx,yptsx,nxy)
178 ALLOCATE(rlatx(no),rlonx(no),xptsx(no),yptsx(no),nxy(no))
191 IF(abs(xpts(n)-fill).GT.tinyreal.AND.abs(ypts(n)-fill).GT.tinyreal)
THEN
192 nxy(n) = grid_in%field_pos(nint(xpts(n)), nint(ypts(n)))
201 IF(iret.EQ.0.AND.iretx.EQ.0)
THEN
202 IF(.not. to_station_points)
THEN
220 IF(ibi(k).EQ.0.OR.li(nxy(n),k))
THEN
224 ELSEIF(mspiral.GT.1)
THEN
227 ixs=int(sign(1.,xpts(n)-i1))
228 jxs=int(sign(1.,ypts(n)-j1))
230 kxs=int(sqrt(4*mx-2.5))
232 SELECT CASE(mod(kxs,4))
234 ix=i1-ixs*(kxs/4-kxt)
238 jx=j1-jxs*(kxs/4-kxt)
240 ix=i1+ixs*(1+kxs/4-kxt)
244 jx=j1+jxs*(kxs/4-kxt)
246 nx = grid_in%field_pos(ix, jx)
261 IF(.NOT.all(lo(1:no,k))) ibo(k)=1
264 select type(grid_out)
265 type is(ip_equid_cylind_grid)
266 CALL polfixs(no,mo,km,rlat,ibo,lo,go)
270 IF(iret.EQ.0) iret=iretx
271 IF(.not. to_station_points) no=0
353 MI,MO,KM,IBI,LI,UI,VI, &
354 NO,RLAT,RLON,CROT,SROT,IBO,LO,UO,VO,IRET)
356 class(ip_grid),
intent(in) :: grid_in, grid_out
357 INTEGER,
INTENT(IN ) :: IPOPT(20)
358 INTEGER,
INTENT(IN ) :: IBI(KM),MI,MO,KM
359 INTEGER,
INTENT(INOUT) :: NO
360 INTEGER,
INTENT( OUT) :: IRET, IBO(KM)
362 LOGICAL*1,
INTENT(IN ) :: LI(MI,KM)
363 LOGICAL*1,
INTENT( OUT) :: LO(MO,KM)
365 REAL,
INTENT(IN ) :: UI(MI,KM),VI(MI,KM)
366 REAL,
INTENT(INOUT) :: CROT(MO),SROT(MO)
367 REAL,
INTENT(INOUT) :: RLAT(MO),RLON(MO)
368 REAL,
INTENT( OUT) :: UO(MO,KM),VO(MO,KM)
370 REAL,
PARAMETER :: FILL=-9999.
372 INTEGER :: I1,J1,IXS,JXS,MX
373 INTEGER :: KXS,KXT,IX,JX,NX
374 INTEGER :: MSPIRAL,N,K,NK,NV
376 LOGICAL :: SAME_GRIDI, SAME_GRIDO
378 REAL :: CX,SX,CM,SM,UROT,VROT
379 REAL :: XPTS(MO),YPTS(MO)
380 REAL :: CROI(MI),SROI(MI)
381 REAL :: XPTI(MI),YPTI(MI),RLOI(MI),RLAI(MI)
383 logical :: to_station_points
385 INTEGER,
SAVE :: NOX=-1,iretx=-1
386 INTEGER,
ALLOCATABLE,
SAVE :: NXY(:)
388 REAL,
ALLOCATABLE,
SAVE :: RLATX(:),RLONX(:),XPTSX(:),YPTSX(:)
389 REAL,
ALLOCATABLE,
SAVE :: CROTX(:),SROTX(:),CXY(:),SXY(:)
390 class(ip_grid),
allocatable,
save :: prev_grid_in, prev_grid_out
394 mspiral=max(ipopt(1),1)
397 if (.not.
allocated(prev_grid_in) .or. .not.
allocated(prev_grid_out))
then
398 allocate(prev_grid_in, source = grid_in)
399 allocate(prev_grid_out, source = grid_out)
404 same_gridi = grid_in == prev_grid_in
405 same_grido = grid_out == prev_grid_out
407 if (.not. same_gridi .or. .not. same_grido)
then
408 deallocate(prev_grid_in)
409 deallocate(prev_grid_out)
411 allocate(prev_grid_in, source = grid_in)
412 allocate(prev_grid_out, source = grid_out)
416 select type(grid_out)
417 type is(ip_station_points_grid)
418 to_station_points = .true.
420 to_station_points = .false.
425 IF(iret.EQ.0.AND.(to_station_points.OR..NOT.same_gridi.OR..NOT.same_grido))
THEN
428 CALL gdswzd(grid_out, 0,mo,fill,xpts,ypts,rlon,rlat,no,crot,srot)
432 CALL gdswzd(grid_in,-1,no,fill,xpts,ypts,rlon,rlat,nv)
433 IF(iret.EQ.0.AND.nv.EQ.0) iret=2
434 CALL gdswzd(grid_in, 0,mi,fill,xpti,ypti,rloi,rlai,nv,croi,sroi)
438 IF(nox.GE.0)
DEALLOCATE(rlatx,rlonx,xptsx,yptsx,crotx,srotx,nxy,cxy,sxy)
439 ALLOCATE(rlatx(no),rlonx(no),xptsx(no),yptsx(no), &
440 crotx(no),srotx(no),nxy(no),cxy(no),sxy(no))
455 IF(abs(xpts(n)-fill).GT.tinyreal.AND.abs(ypts(n)-fill).GT.tinyreal)
THEN
456 nxy(n) = grid_in%field_pos(nint(xpts(n)),nint(ypts(n)))
458 CALL movect(rlai(nxy(n)),rloi(nxy(n)),rlat(n),rlon(n),cm,sm)
459 cxy(n)=cm*croi(nxy(n))+sm*sroi(nxy(n))
460 sxy(n)=sm*croi(nxy(n))-cm*sroi(nxy(n))
470 IF(iret.EQ.0.AND.iretx.EQ.0)
THEN
471 IF(.not. to_station_points)
THEN
494 IF(ibi(k).EQ.0.OR.li(nxy(n),k))
THEN
495 urot=cxy(n)*ui(nxy(n),k)-sxy(n)*vi(nxy(n),k)
496 vrot=sxy(n)*ui(nxy(n),k)+cxy(n)*vi(nxy(n),k)
497 uo(n,k)=crot(n)*urot-srot(n)*vrot
498 vo(n,k)=srot(n)*urot+crot(n)*vrot
501 ELSEIF(mspiral.GT.1)
THEN
504 ixs=int(sign(1.,xpts(n)-i1))
505 jxs=int(sign(1.,ypts(n)-j1))
507 kxs=int(sqrt(4*mx-2.5))
509 SELECT CASE(mod(kxs,4))
511 ix=i1-ixs*(kxs/4-kxt)
515 jx=j1-jxs*(kxs/4-kxt)
517 ix=i1+ixs*(1+kxs/4-kxt)
521 jx=j1+jxs*(kxs/4-kxt)
523 nx = grid_in%field_pos(ix, jx)
526 CALL movect(rlai(nx),rloi(nx),rlat(n),rlon(n),cm,sm)
527 cx=cm*croi(nx)+sm*sroi(nx)
528 sx=sm*croi(nx)-cm*sroi(nx)
529 urot=cx*ui(nx,k)-sx*vi(nx,k)
530 vrot=sx*ui(nx,k)+cx*vi(nx,k)
531 uo(n,k)=crot(n)*urot-srot(n)*vrot
532 vo(n,k)=srot(n)*urot+crot(n)*vrot
543 IF(.NOT.all(lo(1:no,k))) ibo(k)=1
546 select type(grid_out)
547 type is(ip_equid_cylind_grid)
548 CALL polfixv(no,mo,km,rlat,rlon,ibo,lo,uo,vo)
553 IF(iret.EQ.0) iret=iretx
554 IF(.not. to_station_points) no=0
subroutine movect(FLAT, FLON, TLAT, TLON, CROT, SROT)
This subprogram provides the rotation parameters to move a vector along a great circle from one posit...
Driver module for gdswzd routines.
Re-export the individual grids.
Interpolate scalar fields (neighbor).
subroutine interpolate_neighbor_vector(IPOPT, grid_in, grid_out, MI, MO, KM, IBI, LI, UI, VI, NO, RLAT, RLON, CROT, SROT, IBO, LO, UO, VO, IRET)
Interpolate vector fields (neighbor).
subroutine interpolate_neighbor_scalar(IPOPT, grid_in, grid_out, MI, MO, KM, IBI, LI, GI, NO, RLAT, RLON, IBO, LO, GO, IRET)
Interpolate scalar fields (neighbor).
Make multiple pole scalar values consistent.
subroutine, public polfixs(NM, NX, KM, RLAT, IB, LO, GO)
Make multiple pole scalar values consistent.
subroutine, public polfixv(NM, NX, KM, RLAT, RLON, IB, LO, UO, VO)
Make multiple pole vector values consistent,.