/[MITgcm]/MITgcm/pkg/fizhi/fizhi_fillnegs.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/fizhi_fillnegs.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.9 - (show annotations) (download)
Thu Apr 2 20:54:03 2009 UTC (15 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64, checkpoint65, checkpoint62, checkpoint63, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint61z, checkpoint61x, checkpoint61y, HEAD
Changes since 1.8: +33 -33 lines
avoid blank after 72 column (produces warning with some compiler)

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_fillnegs.F,v 1.8 2004/08/05 17:06:40 molod Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 SUBROUTINE QCHECK ( idim1,idim2,jdim1,jdim2,ldim,Nsx,Nsy,im1,im2,
6 . jm1,jm2,bi,bj,dp,qz)
7 C***********************************************************************
8 C Purpose
9 C Check Specific Humidity Field for Negative values
10 C
11 C Argument Description
12 C IDIM1 .... Start Zonal Dimension
13 C IDIM2 .... End Zonal Dimension
14 C JDIM1 .... Start Meridional Dimension
15 C JDIM2 .... End Meridional Dimension
16 C IM1 ...... Start Zonal Span
17 C IM2 ...... End Zonal Span
18 C JM1 ...... Start Meridional Span
19 C JM2 ...... End Meridional Span
20 C LDIM ..... Vertical Dimension
21 C DP ....... Delta Pressure
22 C QZ ........Specific Humidity (g/g)
23 C***********************************************************************
24 implicit none
25
26 integer idim1,idim2,jdim1,jdim2,Ldim,im1,im2,jm1,jm2
27 integer Nsx,Nsy,bi,bj
28 _RL qz(idim1:idim2,jdim1:jdim2,ldim,Nsx,Nsy)
29 _RL dp(idim1:idim2,jdim1:jdim2,ldim,Nsx,Nsy)
30
31 integer i,j,L,LM1
32 _RL ddsig
33
34 c Fill Negative Specific Humidities
35 c ---------------------------------
36 DO L=2,Ldim
37 LM1 = L-1
38 do j=jm1,jm2
39 do i=im1,im2
40 ddsig = dp(i,j,LM1,bi,bj)/dp(i,j,L,bi,bj)
41 if( qz(i,j,LM1,bi,bj).lt.0.0 _d 0) then
42 qz(i,j,L,bi,bj ) = qz(i,j,L,bi,bj) + qz(i,j,LM1,bi,bj)*ddsig
43 qz(i,j,LM1,bi,bj) = 0.0 _d 0
44 endif
45 enddo
46 enddo
47 enddo
48
49 do j=jm1,jm2
50 do i=im1,im2
51 if(qz(i,j,Ldim,bi,bj).lt.0.0 _d 0)qz(i,j,Ldim,bi,bj) = 0.0 _d 0
52 enddo
53 enddo
54
55 return
56 end
57
58 subroutine tracer_fill ( pq,im,jm,lm,dlam,dphi,dp)
59 C***********************************************************************
60 C PURPOSE
61 C Fill negative tracer values using local borrowing
62 C
63 C INPUT
64 C pq ..... Mass-weighted (PI) Tracer
65 C im ..... Zonal Dimension
66 C jm ..... Meridional Dimension
67 C lm ..... Vertical Dimension
68 C dlam ... Zonal Grid Increment
69 C dphi ... Meridional Grid Increment
70 C dp ..... Vertical Grid Increment
71 C
72 C Note:
73 C If no immediate surrounding value is large enough to fill negative
74 C value,
75 C the sum of immediate surrounding positive values is tried.
76 C If sum is not large enough, tracer is simply set to zero.
77 C
78 C***********************************************************************
79
80 implicit none
81
82 c Input Variables
83 c ---------------
84 integer im,jm,lm
85 _RL pq(im,jm,lm),dlam(im),dphi(jm),dp(im,jm,lm)
86
87 c Local Variables
88 c ---------------
89 integer i,j,l,im1,ip1,imax,m
90 _RL lam(im), phi(jm)
91 _RL array(6)
92 _RL pi,a,getcon,undef
93 _RL qmax,qval,sum,fact
94
95 _RL dxu(im,jm)
96 _RL dxv(im,jm)
97 _RL dxp(im,jm)
98 _RL dyv(im,jm)
99 _RL dyp(im,jm)
100
101 _RL d2p(im,jm)
102
103 C *********************************************************
104 C **** Initialization ****
105 C *********************************************************
106
107 pi = 4.0*atan(1.0)
108 a = getcon('EARTH RADIUS')
109
110 c Compute Longitudes
111 c ------------------
112 lam(1) = -pi
113 do i=2,im
114 lam(i) = lam(i-1) + dlam(i-1)
115 enddo
116
117 c Compute Latitudes
118 c -----------------
119 phi(1) = -pi/2.
120 do j=2,jm-1
121 phi(j) = phi(j-1) + dphi(j-1)
122 enddo
123 phi(jm) = pi/2.
124
125 c Compute DXU and DYV
126 c -------------------
127 do j=2,jm-1
128 do i=1,im
129 dxu(i,j) = a*cos(phi(j))*dlam(i)
130 enddo
131 enddo
132
133 do j=2,jm-2
134 do i=1,im
135 dyv(i,j) = a*dphi(j)
136 enddo
137 enddo
138 do i=1,im
139 dyv(i,1) = a*(dphi(1) +0.5*dphi(2) )
140 dyv(i,jm-1) = a*(dphi(jm-1)+0.5*dphi(jm-2))
141 enddo
142
143 c Compute DXP and DXV
144 c -------------------
145 do j=2,jm-1
146 im1 = im
147 do i=1,im
148 dxp(i,j) = ( dxu(i,j)+dxu(im1,j) )*0.5
149 im1 = i
150 enddo
151 enddo
152
153 do j=2,jm-2
154 do i=1,im
155 dxv(i,j) = ( dxp(i,j)+dxp(i,j+1) )*0.5
156 enddo
157 enddo
158
159 c Compute DYP
160 c -----------
161 do j=3,jm-2
162 do i=1,im
163 dyp(i,j) = ( dyv(i,j)+dyv(i,j-1) )*0.5
164 enddo
165 enddo
166 do i=1,im
167 dyp(i,2) = dyv(i,1)
168 dyp(i,jm-1) = dyv(i,jm-1)
169 enddo
170
171 c Compute Area Factor D2P
172 c -----------------------
173 do j=3,jm-2
174 do i=1,im
175 d2p(i,j) = 0.5*( dxv(i,j)+dxv(i,j-1) )*dyp(i,j)
176 enddo
177 enddo
178 do i=1,im
179 d2p(i,2) = dxv(i,2) *dyp(i,2)
180 d2p(i,jm-1) = dxv(i,jm-2)*dyp(i,jm-1)
181 enddo
182
183 undef = getcon('UNDEF')
184
185 C *********************************************************
186 C **** Fill Negative Values ****
187 C *********************************************************
188
189 do l=1,lm
190 do j=2,jm-1
191
192 im1 = im-1
193 i = im
194 do ip1=1,im
195
196 if( pq(i,j,L).lt.0.0 ) then
197
198 qval = pq(i ,j,L)*d2p(i ,j)*dp(i,j,L)
199 array(1) = pq(ip1,j,L)*d2p(ip1,j)*dp(i,j,L)
200 array(2) = pq(im1,j,L)*d2p(im1,j)*dp(i,j,L)
201
202 if( j.eq.jm-1 ) then
203 array(3) = -undef
204 else
205 array(3) = pq(i,j+1,L)*d2p(i,j+1)*dp(i,j,L)
206 endif
207 if( j.eq.2 ) then
208 array(4) = -undef
209 else
210 array(4) = pq(i,j-1,L)*d2p(i,j-1)*dp(i,j,L)
211 endif
212 if( L.eq.1 ) then
213 array(5) = -undef
214 else
215 array(5) = pq(i,j,L-1)*d2p(i,j)*dp(i,j,L)
216 endif
217 if( L.eq.lm ) then
218 array(6) = -undef
219 else
220 array(6) = pq(i,j,L+1)*d2p(i,j)*dp(i,j,L)
221 endif
222
223 call maxval1 (array,6,-qval,qmax,imax)
224
225 if( imax.eq.0 ) then
226 sum = 0.0
227 do m=1,6
228 if( array(m).gt.0.0 ) sum = sum + array(m)
229 enddo
230 if( sum.gt.-qval ) then
231 fact = 1.0 + qval/sum
232 if( array(1).gt.0 ) pq(ip1,j,L) = pq(ip1,j,L) * fact
233 if( array(2).gt.0 ) pq(im1,j,L) = pq(im1,j,L) * fact
234 if( array(3).gt.0 ) pq(i,j+1,L) = pq(i,j+1,L) * fact
235 if( array(4).gt.0 ) pq(i,j-1,L) = pq(i,j-1,L) * fact
236 if( array(5).gt.0 ) pq(i,j,L-1) = pq(i,j,L-1) * fact
237 if( array(6).gt.0 ) pq(i,j,L+1) = pq(i,j,L+1) * fact
238 pq(i,j,L) = 0.0
239 else
240 pq(i,j,L) = 0.0
241 endif
242 else
243 if( imax.eq.1 ) pq(ip1,j,L) = pq(ip1,j,L) +
244 . pq(i,j,L)*d2p(i,j)/d2p(ip1,j)
245 if( imax.eq.2 ) pq(im1,j,L) = pq(im1,j,L) +
246 . pq(i,j,L)*d2p(i,j)/d2p(im1,j)
247 if( imax.eq.3 ) pq(i,j+1,L) = pq(i,j+1,L) +
248 . pq(i,j,L)*d2p(i,j)/d2p(i,j+1)
249 if( imax.eq.4 ) pq(i,j-1,L) = pq(i,j-1,L) +
250 . pq(i,j,L)*d2p(i,j)/d2p(i,j-1)
251 if( imax.eq.5 ) pq(i,j,L-1) = pq(i,j,L-1) +
252 . pq(i,j,L)*dp(i,j,L) /dp(i,j,L-1)
253 if( imax.eq.6 ) pq(i,j,L+1) = pq(i,j,L+1) +
254 . pq(i,j,L)*dp(i,j,L) /dp(i,j,L+1)
255 pq(i,j,L) = 0.0
256 endif
257
258 endif ! End pq<0 Test
259
260 im1 = i
261 i = ip1
262 enddo
263 enddo
264 enddo
265
266 return
267 end
268
269 subroutine maxval1 (q,im,qval,qmax,imax)
270 C***********************************************************************
271 C PURPOSE
272 C Find the location and value of the array element which is greater
273 C than a prescribed value.
274 C
275 C INPUT
276 C q ...... Array Elements
277 C im ..... Dimension of Array q
278 C qval ... Prescribed Value
279 C
280 C OUTPUT
281 C qmax ... Largest Array element which is greater than qval
282 C imax ... Location of Largest Array Element
283 C
284 C Note:
285 C If no array element is larger than qval, then imax = 0
286 C
287 C***********************************************************************
288 implicit none
289 integer im, i, imax
290 _RL q(im), qmax, qval
291 qmax = qval
292 imax = 0
293 do i=1,im
294 if( q(i).gt.qmax ) then
295 qmax = q(i)
296 imax = i
297 endif
298 enddo
299 return
300 end

  ViewVC Help
Powered by ViewVC 1.1.22