1 |
jmc |
1.2 |
! $Header: /u/gcmpack/MITgcm/pkg/atm_phys/lscale_cond_mod.F90,v 1.1 2013/05/08 22:14:14 jmc Exp $ |
2 |
jmc |
1.1 |
! $Name: $ |
3 |
|
|
|
4 |
|
|
module lscale_cond_mod |
5 |
|
|
|
6 |
|
|
!----------------------------------------------------------------------- |
7 |
|
|
!use fms_mod, only: file_exist, error_mesg, open_file, & |
8 |
|
|
! check_nml_error, mpp_pe, FATAL, & |
9 |
|
|
! close_file |
10 |
|
|
use simple_sat_vapor_pres_mod, only: escomp, descomp |
11 |
|
|
|
12 |
|
|
use constants_mod, only: HLv,HLs,Cp_air,Grav,rdgas,rvgas |
13 |
|
|
use gcm_params_mod, only: gcm_LEN_MBUF, gcm_SQZ_R, gcm_stdMsgUnit |
14 |
|
|
|
15 |
|
|
implicit none |
16 |
|
|
private |
17 |
|
|
!----------------------------------------------------------------------- |
18 |
|
|
! ---- public interfaces ---- |
19 |
|
|
|
20 |
|
|
public lscale_cond, lscale_cond_init |
21 |
|
|
|
22 |
|
|
!----------------------------------------------------------------------- |
23 |
|
|
! ---- version number ---- |
24 |
|
|
|
25 |
jmc |
1.2 |
character(len=128) :: version = '$Id: lscale_cond_mod.F90,v 1.1 2013/05/08 22:14:14 jmc Exp $' |
26 |
jmc |
1.1 |
character(len=128) :: tag = '$Name: $' |
27 |
|
|
|
28 |
|
|
!----------------------------------------------------------------------- |
29 |
|
|
! ---- local/private data ---- |
30 |
|
|
|
31 |
|
|
real, parameter :: d622 = rdgas/rvgas |
32 |
|
|
real, parameter :: d378 = 1.-d622 |
33 |
|
|
|
34 |
|
|
logical :: do_init=.true. |
35 |
|
|
|
36 |
|
|
!----------------------------------------------------------------------- |
37 |
|
|
! --- namelist ---- |
38 |
|
|
|
39 |
|
|
real :: hc=1.00 |
40 |
|
|
logical :: do_evap=.false. |
41 |
|
|
|
42 |
|
|
namelist /lscale_cond_nml/ hc, do_evap |
43 |
|
|
|
44 |
|
|
!----------------------------------------------------------------------- |
45 |
|
|
! description of namelist variables |
46 |
|
|
! |
47 |
|
|
! hc = relative humidity at which large scale condensation |
48 |
|
|
! occurs, where 0 <= hc <= 1 (default: hc=1.) |
49 |
|
|
! |
50 |
|
|
! do_evap = flag for the re-evaporation of moisture in |
51 |
|
|
! sub-saturated layers below, if do_evap=.true. then |
52 |
|
|
! re-evaporation is performed (default: do_evap=.false.) |
53 |
|
|
! |
54 |
|
|
!----------------------------------------------------------------------- |
55 |
|
|
|
56 |
|
|
contains |
57 |
|
|
|
58 |
|
|
!####################################################################### |
59 |
|
|
|
60 |
|
|
! subroutine lscale_cond (tin, qin, pfull, phalf, coldT, & |
61 |
|
|
! rain, snow, tdel, qdel, mask, conv) |
62 |
|
|
subroutine lscale_cond (tin, qin, pfull, phalf, coldT, & |
63 |
|
|
rain, snow, tdel, qdel, qsat, & |
64 |
|
|
myThid, mask, conv ) |
65 |
|
|
|
66 |
|
|
!----------------------------------------------------------------------- |
67 |
|
|
! |
68 |
|
|
! large scale condensation |
69 |
|
|
! |
70 |
|
|
!----------------------------------------------------------------------- |
71 |
|
|
! |
72 |
|
|
! input: tin temperature at full model levels |
73 |
|
|
! qin specific humidity of water vapor at full |
74 |
|
|
! model levels |
75 |
|
|
! pfull pressure at full model levels |
76 |
|
|
! phalf pressure at half (interface) model levels |
77 |
|
|
! coldT should precipitation be snow at this point? |
78 |
|
|
! optional: |
79 |
|
|
! mask optional mask (0 or 1.) |
80 |
|
|
! conv logical flag; if true then no large-scale |
81 |
|
|
! adjustment is performed at that grid-point or |
82 |
|
|
! model level |
83 |
|
|
! |
84 |
|
|
! output: rain liquid precipitation (kg/m2) |
85 |
|
|
! snow frozen precipitation (kg/m2) |
86 |
|
|
! tdel temperature tendency at full model levels |
87 |
|
|
! qdel specific humidity tendency (of water vapor) at |
88 |
|
|
! full model levels |
89 |
|
|
! qsat saturated specific humidity |
90 |
|
|
! |
91 |
|
|
!----------------------------------------------------------------------- |
92 |
|
|
!--------------------- interface arguments ----------------------------- |
93 |
|
|
|
94 |
|
|
real , intent(in) , dimension(:,:,:) :: tin, qin, pfull, phalf |
95 |
|
|
logical , intent(in) , dimension(:,:):: coldT |
96 |
|
|
real , intent(out), dimension(:,:) :: rain,snow |
97 |
|
|
real , intent(out), dimension(:,:,:) :: tdel, qdel, qsat |
98 |
|
|
integer, intent(in) :: myThid |
99 |
|
|
real , intent(in) , dimension(:,:,:), optional :: mask |
100 |
|
|
logical, intent(in) , dimension(:,:,:), optional :: conv |
101 |
|
|
!----------------------------------------------------------------------- |
102 |
|
|
!---------------------- local data ------------------------------------- |
103 |
|
|
|
104 |
|
|
logical,dimension(size(tin,1),size(tin,2),size(tin,3)) :: do_adjust |
105 |
|
|
real,dimension(size(tin,1),size(tin,2),size(tin,3)) :: & |
106 |
|
|
esat, desat, dqsat, pmes, pmass |
107 |
|
|
real,dimension(size(tin,1),size(tin,2)) :: hlcp, precip |
108 |
|
|
integer k, kx |
109 |
|
|
!----------------------------------------------------------------------- |
110 |
|
|
! computation of precipitation by condensation processes |
111 |
|
|
!----------------------------------------------------------------------- |
112 |
|
|
|
113 |
|
|
! if (do_init) call error_mesg ('lscale_cond', & |
114 |
|
|
! 'lscale_cond_init has not been called.', FATAL) |
115 |
|
|
|
116 |
|
|
kx=size(tin,3) |
117 |
|
|
|
118 |
|
|
!----- compute proper latent heat -------------------------------------- |
119 |
|
|
WHERE (coldT) |
120 |
|
|
hlcp = HLs/Cp_air |
121 |
|
|
ELSEWHERE |
122 |
|
|
hlcp = HLv/Cp_air |
123 |
|
|
END WHERE |
124 |
|
|
|
125 |
|
|
!----- saturation vapor pressure (esat) & specific humidity (qsat) ----- |
126 |
|
|
|
127 |
|
|
call escomp (tin,esat) |
128 |
|
|
call descomp (tin,desat) |
129 |
|
|
|
130 |
|
|
esat(:,:,:)=esat(:,:,:)*hc |
131 |
|
|
|
132 |
|
|
where (pfull(:,:,:) > d378*esat(:,:,:)) |
133 |
|
|
pmes(:,:,:)=1.0/(pfull(:,:,:)-d378*esat(:,:,:)) |
134 |
|
|
qsat(:,:,:)=d622*esat(:,:,:)*pmes(:,:,:) |
135 |
|
|
qsat(:,:,:)=max(0.0,qsat(:,:,:)) |
136 |
|
|
dqsat(:,:,:)=d622*pfull(:,:,:)*desat(:,:,:)*pmes(:,:,:)*pmes(:,:,:) |
137 |
|
|
elsewhere |
138 |
|
|
pmes(:,:,:)=0.0 |
139 |
|
|
qsat(:,:,:)=0.0 |
140 |
|
|
dqsat(:,:,:)=0.0 |
141 |
|
|
endwhere |
142 |
|
|
|
143 |
|
|
!--------- do adjustment where greater than saturated value ------------ |
144 |
|
|
|
145 |
|
|
if (present(conv)) then |
146 |
|
|
!!!! do_adjust(:,:,:)=(.not.conv(:,:,:) .and. qin(:,:,:) > qsat(:,:,:)) |
147 |
|
|
do_adjust(:,:,:)=(.not.conv(:,:,:) .and. & |
148 |
|
|
(qin(:,:,:) - qsat(:,:,:))*qsat(:,:,:) > 0.0) |
149 |
|
|
else |
150 |
|
|
!!!! do_adjust(:,:,:)=(qin(:,:,:) > qsat(:,:,:)) |
151 |
|
|
do_adjust(:,:,:)=( (qin(:,:,:) - qsat(:,:,:))*qsat(:,:,:) > 0.0) |
152 |
|
|
endif |
153 |
|
|
|
154 |
|
|
if (present(mask)) then |
155 |
|
|
do_adjust(:,:,:)=do_adjust(:,:,:) .and. (mask(:,:,:) > 0.5) |
156 |
|
|
end if |
157 |
|
|
|
158 |
|
|
!----------- compute adjustments to temp and spec humidity ------------- |
159 |
|
|
do k = 1,kx |
160 |
|
|
where (do_adjust(:,:,k)) |
161 |
|
|
qdel(:,:,k)=(qsat(:,:,k)-qin(:,:,k))/(1.0+hlcp(:,:)*dqsat(:,:,k)) |
162 |
|
|
tdel(:,:,k)=-hlcp(:,:)*qdel(:,:,k) |
163 |
|
|
elsewhere |
164 |
|
|
qdel(:,:,k)=0.0 |
165 |
|
|
tdel(:,:,k)=0.0 |
166 |
|
|
endwhere |
167 |
|
|
end do |
168 |
|
|
|
169 |
|
|
!------------ pressure mass of each layer ------------------------------ |
170 |
|
|
|
171 |
|
|
do k=1,kx |
172 |
|
|
pmass(:,:,k)=(phalf(:,:,k+1)-phalf(:,:,k))/Grav |
173 |
|
|
enddo |
174 |
|
|
|
175 |
|
|
!------------ re-evaporation of precipitation in dry layer below ------- |
176 |
|
|
|
177 |
|
|
if (do_evap) then |
178 |
|
|
if (present(mask)) then |
179 |
|
|
! call precip_evap (pmass,tin,qin,qsat,dqsat,hlcp,tdel,qdel,mask) |
180 |
|
|
call precip_evap (pmass,tin,qin,qsat,dqsat,hlcp,tdel,qdel, & |
181 |
|
|
myThid, mask ) |
182 |
|
|
else |
183 |
|
|
call precip_evap (pmass,tin,qin,qsat,dqsat,hlcp,tdel,qdel, & |
184 |
|
|
myThid ) |
185 |
|
|
endif |
186 |
|
|
endif |
187 |
|
|
|
188 |
|
|
!------------ integrate precip ----------------------------------------- |
189 |
|
|
|
190 |
|
|
precip(:,:)=0.0 |
191 |
|
|
do k=1,kx |
192 |
|
|
precip(:,:)=precip(:,:)-pmass(:,:,k)*qdel(:,:,k) |
193 |
|
|
enddo |
194 |
|
|
precip(:,:)=max(precip(:,:),0.0) |
195 |
|
|
|
196 |
|
|
!assign precip to snow or rain |
197 |
|
|
WHERE (coldT) |
198 |
|
|
snow = precip |
199 |
|
|
rain = 0. |
200 |
|
|
ELSEWHERE |
201 |
|
|
rain = precip |
202 |
|
|
snow = 0. |
203 |
|
|
END WHERE |
204 |
|
|
|
205 |
|
|
!----------------------------------------------------------------------- |
206 |
|
|
|
207 |
|
|
end subroutine lscale_cond |
208 |
|
|
|
209 |
|
|
!####################################################################### |
210 |
|
|
|
211 |
|
|
! subroutine precip_evap (pmass, tin, qin, qsat, dqsat, hlcp, & |
212 |
|
|
! tdel, qdel, mask) |
213 |
|
|
subroutine precip_evap (pmass, tin, qin, qsat, dqsat, hlcp, & |
214 |
|
|
tdel, qdel, myThid, mask ) |
215 |
|
|
|
216 |
|
|
!----------------------------------------------------------------------- |
217 |
|
|
! performs re-evaporation of falling precipitation |
218 |
|
|
!----------------------------------------------------------------------- |
219 |
|
|
real, intent(in), dimension(:,:,:) :: pmass, tin, qin, qsat, dqsat |
220 |
|
|
real, intent(in), dimension(:,:) :: hlcp |
221 |
|
|
real, intent(inout), dimension(:,:,:) :: tdel, qdel |
222 |
|
|
integer, intent(in) :: myThid |
223 |
|
|
real, intent(in), dimension(:,:,:), optional :: mask |
224 |
|
|
!----------------------------------------------------------------------- |
225 |
|
|
real, dimension(size(tin,1),size(tin,2)) :: exq, def |
226 |
|
|
|
227 |
|
|
integer k |
228 |
|
|
!----------------------------------------------------------------------- |
229 |
|
|
exq(:,:)=0.0 |
230 |
|
|
|
231 |
|
|
do k=1,size(tin,3) |
232 |
|
|
|
233 |
|
|
where (qdel(:,:,k) < 0.0) exq(:,:) = exq(:,:) - & |
234 |
|
|
qdel(:,:,k)*pmass(:,:,k) |
235 |
|
|
|
236 |
|
|
if (present(mask)) exq(:,:) = exq(:,:)*mask(:,:,k) |
237 |
|
|
|
238 |
|
|
! ---- evaporate precip where needed ------ |
239 |
|
|
|
240 |
|
|
where ( (qdel(:,:,k) >= 0.0) .and. (exq(:,:) > 0.0) ) |
241 |
|
|
exq(:,:) = exq(:,:) / pmass(:,:,k) |
242 |
|
|
def(:,:) = (qsat(:,:,k)-qin(:,:,k))/(1.+hlcp(:,:)*dqsat(:,:,k)) |
243 |
|
|
def(:,:) = min(max(def(:,:),0.0),exq(:,:)) |
244 |
|
|
qdel(:,:,k) = qdel(:,:,k) + def(:,:) |
245 |
|
|
tdel(:,:,k) = tdel(:,:,k) - def(:,:)*hlcp(:,:) |
246 |
|
|
exq(:,:) = (exq(:,:)-def(:,:))*pmass(:,:,k) |
247 |
|
|
endwhere |
248 |
|
|
|
249 |
|
|
enddo |
250 |
|
|
|
251 |
|
|
!----------------------------------------------------------------------- |
252 |
|
|
|
253 |
|
|
end subroutine precip_evap |
254 |
|
|
|
255 |
|
|
!####################################################################### |
256 |
|
|
|
257 |
|
|
subroutine lscale_cond_init (myThid) |
258 |
|
|
|
259 |
|
|
!----------------------------------------------------------------------- |
260 |
|
|
! |
261 |
|
|
! initialization for large scale condensation |
262 |
|
|
! |
263 |
|
|
!----------------------------------------------------------------------- |
264 |
|
|
|
265 |
|
|
integer, intent(in) :: myThid |
266 |
|
|
|
267 |
|
|
!----------------------------------------------------------------------- |
268 |
|
|
|
269 |
|
|
! integer unit,io,ierr |
270 |
|
|
integer :: iUnit |
271 |
|
|
CHARACTER*(gcm_LEN_MBUF) :: msgBuf |
272 |
|
|
|
273 |
|
|
!----------- read namelist --------------------------------------------- |
274 |
|
|
|
275 |
|
|
! _BARRIER |
276 |
|
|
! _BEGIN_MASTER(myThid) |
277 |
|
|
CALL BARRIER(myThid) |
278 |
|
|
IF ( myThid.EQ.1 ) THEN |
279 |
|
|
|
280 |
|
|
WRITE(msgBuf,'(A)') 'LSCALE_COND_INIT: opening data.atm_gray' |
281 |
|
|
CALL PRINT_MESSAGE( msgBuf, gcm_stdMsgUnit, gcm_SQZ_R, myThid ) |
282 |
|
|
CALL OPEN_COPY_DATA_FILE( & |
283 |
|
|
'data.atm_gray', 'LSCALE_COND_INIT', & |
284 |
|
|
iUnit, & |
285 |
|
|
myThid ) |
286 |
|
|
! Read parameters from open data file |
287 |
|
|
READ(UNIT=iUnit,NML=lscale_cond_nml) |
288 |
|
|
WRITE(msgBuf,'(A)') & |
289 |
|
|
'LSCALE_COND_INIT: finished reading data.atm_gray' |
290 |
|
|
CALL PRINT_MESSAGE( msgBuf, gcm_stdMsgUnit, gcm_SQZ_R, myThid ) |
291 |
|
|
! Close the open data file |
292 |
jmc |
1.2 |
#ifdef SINGLE_DISK_IO |
293 |
jmc |
1.1 |
CLOSE(iUnit) |
294 |
jmc |
1.2 |
#else |
295 |
|
|
CLOSE(iUnit,STATUS='DELETE') |
296 |
|
|
#endif /* SINGLE_DISK_IO */ |
297 |
jmc |
1.1 |
|
298 |
|
|
! if (file_exist('input.nml')) then |
299 |
|
|
! unit = open_file (file='input.nml', action='read') |
300 |
|
|
! ierr=1; do while (ierr /= 0) |
301 |
|
|
! read (unit, nml=lscale_cond_nml, iostat=io, end=10) |
302 |
|
|
! ierr = check_nml_error (io,'lscale_cond_nml') |
303 |
|
|
! enddo |
304 |
|
|
! 10 call close_file (unit) |
305 |
|
|
! endif |
306 |
|
|
|
307 |
|
|
!---------- output namelist -------------------------------------------- |
308 |
|
|
|
309 |
|
|
! unit = open_file (file='logfile.out', action='append') |
310 |
|
|
! if ( mpp_pe() == 0 ) then |
311 |
|
|
! write (unit,'(/,80("="),/(a))') trim(version), trim(tag) |
312 |
|
|
! write (unit,nml=lscale_cond_nml) |
313 |
|
|
! endif |
314 |
|
|
! call close_file (unit) |
315 |
|
|
|
316 |
|
|
do_init=.false. |
317 |
|
|
|
318 |
|
|
ENDIF |
319 |
|
|
CALL BARRIER(myThid) |
320 |
|
|
|
321 |
|
|
end subroutine lscale_cond_init |
322 |
|
|
|
323 |
|
|
!####################################################################### |
324 |
|
|
|
325 |
|
|
end module lscale_cond_mod |
326 |
|
|
|