NCEPLIBS-ip  4.3.0
polfix_mod.F90
Go to the documentation of this file.
1 
4 
7 module polfix_mod
8  implicit none
9 
10  private
11  public :: polfixs, polfixv
12 
13 contains
14 
29  SUBROUTINE polfixs(NM,NX,KM,RLAT,IB,LO,GO)
30  IMPLICIT NONE
31  !
32  INTEGER, INTENT(IN ) :: nm, nx, km
33  INTEGER, INTENT(IN ) :: ib(km)
34  !
35  LOGICAL*1, INTENT(INOUT) :: lo(nx,km)
36  !
37  REAL, INTENT(IN ) :: rlat(nm)
38  REAL, INTENT(INOUT) :: go(nx,km)
39  !
40  REAL, PARAMETER :: rlatnp=89.9995
41  REAL, PARAMETER :: rlatsp=-rlatnp
42  !
43  INTEGER :: k, n
44  !
45  REAL :: wnp, gnp, tnp, wsp, gsp, tsp
46  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47  DO k=1,km
48  wnp=0.
49  gnp=0.
50  tnp=0.
51  wsp=0.
52  gsp=0.
53  tsp=0.
54  ! AVERAGE MULTIPLE POLE VALUES
55  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:WNP,GNP,TNP,WSP,GSP,TSP) SCHEDULE(STATIC)
56  DO n=1,nm
57  IF(rlat(n).GE.rlatnp) THEN
58  wnp=wnp+1
59  IF(ib(k).EQ.0.OR.lo(n,k)) THEN
60  gnp=gnp+go(n,k)
61  tnp=tnp+1
62  ENDIF
63  ELSEIF(rlat(n).LE.rlatsp) THEN
64  wsp=wsp+1
65  IF(ib(k).EQ.0.OR.lo(n,k)) THEN
66  gsp=gsp+go(n,k)
67  tsp=tsp+1
68  ENDIF
69  ENDIF
70  ENDDO
71  !$OMP END PARALLEL DO
72  ! DISTRIBUTE AVERAGE VALUES BACK TO MULTIPLE POLES
73  IF(wnp.GT.1) THEN
74  IF(tnp.GE.wnp/2) THEN
75  gnp=gnp/tnp
76  ELSE
77  gnp=0.
78  ENDIF
79  !$OMP PARALLEL DO PRIVATE(N) SCHEDULE(STATIC)
80  DO n=1,nm
81  IF(rlat(n).GE.rlatnp) THEN
82  IF(ib(k).NE.0) lo(n,k)=tnp.GE.wnp/2
83  go(n,k)=gnp
84  ENDIF
85  ENDDO
86  !$OMP END PARALLEL DO
87  ENDIF
88  IF(wsp.GT.1) THEN
89  IF(tsp.GE.wsp/2) THEN
90  gsp=gsp/tsp
91  ELSE
92  gsp=0.
93  ENDIF
94  !$OMP PARALLEL DO PRIVATE(N) SCHEDULE(STATIC)
95  DO n=1,nm
96  IF(rlat(n).LE.rlatsp) THEN
97  IF(ib(k).NE.0) lo(n,k)=tsp.GE.wsp/2
98  go(n,k)=gsp
99  ENDIF
100  ENDDO
101  !$OMP END PARALLEL DO
102  ENDIF
103  ENDDO
104  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105  END SUBROUTINE polfixs
106 
124  SUBROUTINE polfixv(NM,NX,KM,RLAT,RLON,IB,LO,UO,VO)
125  IMPLICIT NONE
126  !
127  INTEGER, INTENT(IN ) :: ib(km), nm, nx, km
128  !
129  LOGICAL*1, INTENT(INOUT) :: lo(nx,km)
130  !
131  REAL, INTENT(IN ) :: rlat(nm), rlon(nm)
132  REAL, INTENT(INOUT) :: uo(nx,km), vo(nx,km)
133  !
134  REAL, PARAMETER :: rlatnp=89.9995
135  REAL, PARAMETER :: rlatsp=-rlatnp
136  REAL, PARAMETER :: pi=3.14159265358979
137  REAL, PARAMETER :: dpr=180./pi
138  !
139  INTEGER :: k, n
140  !
141  REAL :: clon(nm),slon(nm)
142  REAL :: tnp, unp, vnp, wnp
143  REAL :: tsp, usp, vsp, wsp
144  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145  !$OMP PARALLEL DO PRIVATE(N) SCHEDULE(STATIC)
146  DO n=1,nm
147  clon(n)=cos(rlon(n)/dpr)
148  slon(n)=sin(rlon(n)/dpr)
149  ENDDO
150  !$OMP END PARALLEL DO
151  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152  DO k=1,km
153  wnp=0.
154  unp=0.
155  vnp=0.
156  tnp=0.
157  wsp=0.
158  usp=0.
159  vsp=0.
160  tsp=0.
161  ! AVERAGE MULTIPLE POLE VALUES
162  !$OMP PARALLEL DO PRIVATE(N) REDUCTION(+:WNP,UNP,VNP,TNP,WSP,USP,VSP,TSP) SCHEDULE(STATIC)
163  DO n=1,nm
164  IF(rlat(n).GE.rlatnp) THEN
165  wnp=wnp+1
166  IF(ib(k).EQ.0.OR.lo(n,k)) THEN
167  unp=unp+(clon(n)*uo(n,k)-slon(n)*vo(n,k))
168  vnp=vnp+(slon(n)*uo(n,k)+clon(n)*vo(n,k))
169  tnp=tnp+1
170  ENDIF
171  ELSEIF(rlat(n).LE.rlatsp) THEN
172  wsp=wsp+1
173  IF(ib(k).EQ.0.OR.lo(n,k)) THEN
174  usp=usp+(clon(n)*uo(n,k)+slon(n)*vo(n,k))
175  vsp=vsp+(-slon(n)*uo(n,k)+clon(n)*vo(n,k))
176  tsp=tsp+1
177  ENDIF
178  ENDIF
179  ENDDO
180  !$OMP END PARALLEL DO
181  ! DISTRIBUTE AVERAGE VALUES BACK TO MULTIPLE POLES
182  IF(wnp.GT.1) THEN
183  IF(tnp.GE.wnp/2) THEN
184  unp=unp/tnp
185  vnp=vnp/tnp
186  ELSE
187  unp=0.
188  vnp=0.
189  ENDIF
190  !$OMP PARALLEL DO PRIVATE(N) SCHEDULE(STATIC)
191  DO n=1,nm
192  IF(rlat(n).GE.rlatnp) THEN
193  IF(ib(k).NE.0) lo(n,k)=tnp.GE.wnp/2
194  uo(n,k)=clon(n)*unp+slon(n)*vnp
195  vo(n,k)=-slon(n)*unp+clon(n)*vnp
196  ENDIF
197  ENDDO
198  !$OMP END PARALLEL DO
199  ENDIF
200  IF(wsp.GT.1) THEN
201  IF(tsp.GE.wsp/2) THEN
202  usp=usp/wsp
203  vsp=vsp/wsp
204  ELSE
205  usp=0.
206  vsp=0.
207  ENDIF
208  !$OMP PARALLEL DO PRIVATE(N) SCHEDULE(STATIC)
209  DO n=1,nm
210  IF(rlat(n).LE.rlatsp) THEN
211  IF(ib(k).NE.0) lo(n,k)=tsp.GE.wsp/2
212  uo(n,k)=clon(n)*usp-slon(n)*vsp
213  vo(n,k)=slon(n)*usp+clon(n)*vsp
214  ENDIF
215  ENDDO
216  !$OMP END PARALLEL DO
217  ENDIF
218  ENDDO
219  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220  END SUBROUTINE polfixv
221 end module polfix_mod
Make multiple pole scalar values consistent.
Definition: polfix_mod.F90:7
subroutine, public polfixs(NM, NX, KM, RLAT, IB, LO, GO)
Make multiple pole scalar values consistent.
Definition: polfix_mod.F90:30
subroutine, public polfixv(NM, NX, KM, RLAT, RLON, IB, LO, UO, VO)
Make multiple pole vector values consistent,.
Definition: polfix_mod.F90:125