39 use mpp_domains_mod
, only: mpp_update_domains, domain2d
40 use platform_mod
, only: kind_phys => r8_kind
51 subroutine fillz(im, km, nq, q, dp)
52 integer,
intent(in):: im
53 integer,
intent(in):: km
54 integer,
intent(in):: nq
55 real ,
intent(in):: dp(im,km)
56 real ,
intent(inout) :: q(im,km,nq)
61 real qup, qly, dup, dq, sum0, sum1, fac
69 if( q(i,k,ic) < 0. )
then 70 q(i,k1,ic) = q(i,k1,ic) + q(i,k,ic)*dp(i,k)/dp(i,k1)
79 if( q(i,k,ic) < 0. )
then 80 q(i,k1,ic) = q(i,k1,ic) + q(i,k,ic)*dp(i,k)/dp(i,k1)
88 if( q(i,1,ic) < 0. )
then 89 q(i,2,ic) = q(i,2,ic) + q(i,1,ic)*dp(i,1)/dp(i,2)
98 if( q(i,k,ic) < 0. )
then 100 if ( q(i,k-1,ic) > 0. )
then 102 dq = min( q(i,k-1,ic)*dp(i,k-1), -q(i,k,ic)*dp(i,k) )
103 q(i,k-1,ic) = q(i,k-1,ic) - dq/dp(i,k-1)
104 q(i,k ,ic) = q(i,k ,ic) + dq/dp(i,k )
106 if ( q(i,k,ic)<0.0 .and. q(i,k+1,ic)>0. )
then 108 dq = min( q(i,k+1,ic)*dp(i,k+1), -q(i,k,ic)*dp(i,k) )
109 q(i,k+1,ic) = q(i,k+1,ic) - dq/dp(i,k+1)
110 q(i,k ,ic) = q(i,k ,ic) + dq/dp(i,k )
119 if( q(i,k,ic)<0. .and. q(i,k-1,ic)>0.)
then 122 qup = q(i,k-1,ic)*dp(i,k-1)
123 qly = -q(i,k ,ic)*dp(i,k )
125 q(i,k-1,ic) = q(i,k-1,ic) - dup/dp(i,k-1)
126 q(i,k, ic) = q(i,k, ic) + dup/dp(i,k )
136 dm(k) = q(i,k,ic)*dp(i,k)
140 if ( sum0 > 0. )
then 143 sum1 = sum1 + max(0., dm(k))
147 q(i,k,ic) = max(0., fac*dm(k)/dp(i,k))
162 subroutine fill_gfs(im, km, pe2, q, q_min)
164 integer,
intent(in):: im, km
165 real(kind=kind_phys),
intent(in):: pe2(im,km+1)
166 real(kind=kind_phys),
intent(in):: q_min
167 real(kind=kind_phys),
intent(inout):: q(im,km)
169 real(kind=kind_phys) :: dp(im,km)
174 dp(i,k) = pe2(i,k) - pe2(i,k+1)
182 if ( q(i,k)<q_min )
then 184 q(i,k1) = q(i,k1) + (q(i,k)-q_min)*dp(i,k)/dp(i,k1)
194 if ( q(i,k)<0.0 )
then 196 q(i,k1) = q(i,k1) + q(i,k)*dp(i,k)/dp(i,k1)
206 subroutine fill2d(is, ie, js, je, ng, km, q, delp, area, domain, bounded_domain, npx, npy)
208 type(domain2d),
intent(INOUT) :: domain
209 integer,
intent(in):: is, ie, js, je, ng, km, npx, npy
210 logical,
intent(IN):: bounded_domain
211 real,
intent(in):: area(is-ng:ie+ng, js-ng:je+ng)
212 real,
intent(in):: delp(is-ng:ie+ng, js-ng:je+ng, km)
213 real,
intent(inout):: q(is-ng:ie+ng, js-ng:je+ng, km)
215 real,
dimension(is-ng:ie+ng, js-ng:je+ng,km):: qt
216 real,
dimension(is:ie+1, js:je):: fx
217 real,
dimension(is:ie, js:je+1):: fy
218 real,
parameter:: dif = 0.25
220 integer :: is1, ie1, js1, je1
222 if (bounded_domain)
then 228 if (ie == npx-1)
then 238 if (je == npy-1)
then 254 qt(i,j,k) = q(i,j,k)*delp(i,j,k)*area(i,j)
258 call mpp_update_domains(qt, domain, whalo=1, ehalo=1, shalo=1, nhalo=1)
265 if( qt(i-1,j,k)*qt(i,j,k)<0. ) fx(i,j) = qt(i-1,j,k) - qt(i,j,k)
271 if( qt(i,j-1,k)*qt(i,j,k)<0. ) fy(i,j) = qt(i,j-1,k) - qt(i,j,k)
276 q(i,j,k) = q(i,j,k)+dif*(fx(i,j)-fx(i+1,j)+fy(i,j)-fy(i,j+1))/(delp(i,j,k)*area(i,j))
subroutine, public fillz(im, km, nq, q, dp)
The subroutine 'fillz' is for mass-conservative filling of nonphysical negative values in the tracers...
subroutine, public fill_gfs(im, km, pe2, q, q_min)
The subroutine 'fill_gfs' is for mass-conservative filling of nonphysical negative values in the trac...
subroutine, public fill2d(is, ie, js, je, ng, km, q, delp, area, domain, bounded_domain, npx, npy)
The subroutine 'fill2D' fills in nonphysical negative values in a single scalar field using a two-dim...