CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
gfdl_cloud_microphys_v3.F90
1
5
6 use gfdl_cloud_microphys_v3_mod, only: gfdl_cloud_microphys_v3_mod_init, &
7 gfdl_cloud_microphys_v3_mod_driver, &
8 gfdl_cloud_microphys_v3_mod_end, &
9 rad_ref, cld_eff_rad
10
11 implicit none
12
13 private
14
15 public gfdl_cloud_microphys_v3_run, gfdl_cloud_microphys_v3_init, gfdl_cloud_microphys_v3_finalize
16
17 logical :: is_initialized = .false.
18
19contains
20
21! -----------------------------------------------------------------------
22! CCPP entry points for gfdl cloud microphysics
23! -----------------------------------------------------------------------
24
31
32 subroutine gfdl_cloud_microphys_v3_init (me, master, nlunit, input_nml_file, logunit, &
33 fn_nml, imp_physics, imp_physics_gfdl, do_shoc, &
34 hydrostatic, errmsg, errflg)
35
36 implicit none
37
38 integer, intent (in) :: me
39 integer, intent (in) :: master
40 integer, intent (in) :: nlunit
41 integer, intent (in) :: logunit
42 character(len=*), intent (in) :: fn_nml
43 character(len=*), intent (in) :: input_nml_file(:)
44 integer, intent( in) :: imp_physics
45 integer, intent( in) :: imp_physics_gfdl
46 logical, intent( in) :: do_shoc
47 logical, intent( in) :: hydrostatic
48 character(len=*), intent(out) :: errmsg
49 integer, intent(out) :: errflg
50
51 ! Initialize CCPP error handling variables
52 errmsg = ''
53 errflg = 0
54
55 if (is_initialized) return
56
57 if (imp_physics/=imp_physics_gfdl) then
58 write(errmsg,'(*(a))') 'Namelist option for microphysics does not match choice in suite definition file'
59 errflg = 1
60 return
61 end if
62
63 if (do_shoc) then
64 write(errmsg,'(*(a))') 'SHOC is not currently compatible with GFDL MP v3'
65 errflg = 1
66 return
67 endif
68
69 call gfdl_cloud_microphys_v3_mod_init(me, master, nlunit, input_nml_file, logunit, fn_nml, hydrostatic, errmsg, errflg)
70
71 is_initialized = .true.
72
73 end subroutine gfdl_cloud_microphys_v3_init
74
75
76! =======================================================================
83 subroutine gfdl_cloud_microphys_v3_finalize(errmsg, errflg)
84
85 implicit none
86
87 character(len=*), intent(out) :: errmsg
88 integer, intent(out) :: errflg
89
90 ! Initialize CCPP error handling variables
91 errmsg = ''
92 errflg = 0
93
94 if (.not.is_initialized) return
95
96 call gfdl_cloud_microphys_v3_mod_end()
97
98 is_initialized = .false.
99
100 end subroutine gfdl_cloud_microphys_v3_finalize
101
118 subroutine gfdl_cloud_microphys_v3_run(fast_mp_consv, &
119 levs, im, rainmin, con_g, con_fvirt, con_rd, con_eps, garea, slmsk, snowd, &
120 gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, gq0_ntsw, gq0_ntgl, gq0_ntclamt, aerfld, &
121 gt0, gu0, gv0, vvl, prsl, phii, del, &
122 rain0, ice0, snow0, graupel0, prcp0, sr, oro, &
123 dtp, hydrostatic, lradar, refl_10cm, &
124 reset, effr_in, rew, rei, rer, res, reg, &
125 cplchm, pfi_lsan, pfl_lsan, con_one, con_p001, con_secinday, errmsg, errflg)
126
127 use machine, only: kind_phys, kind_dyn, kind_dbl_prec
128
129 implicit none
130
131 ! interface variables
132 integer, intent(in ) :: levs, im
133 real(kind=kind_phys), intent(in ) :: con_g, con_fvirt, con_rd, con_eps, rainmin
134 real(kind=kind_phys), intent(in ) :: con_one, con_p001, con_secinday
135 real(kind=kind_phys), intent(in ), dimension(:) :: garea, slmsk, snowd, oro
136 real(kind=kind_phys), intent(inout), dimension(:,:) :: gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, &
137 gq0_ntsw, gq0_ntgl, gq0_ntclamt
138 real(kind_phys), intent(in ), dimension(:,:,:) :: aerfld
139 real(kind=kind_phys), intent(inout), dimension(:,:) :: gt0, gu0, gv0
140 real(kind=kind_phys), intent(in ), dimension(:,:) :: vvl, prsl, del
141 real(kind=kind_phys), intent(in ), dimension(:,:) :: phii
142
143 ! rain/snow/ice/graupel/precip amounts, fraction of frozen precip
144 !real(kind_phys), dimension(:) :: water0
145 real(kind_phys), intent(out ), dimension(:), optional :: rain0
146 real(kind_phys), intent(out ), dimension(:), optional :: snow0
147 real(kind_phys), intent(out ), dimension(:), optional :: ice0
148 real(kind_phys), intent(out ), dimension(:), optional :: graupel0
149 real(kind_phys), intent(out ), dimension(:) :: prcp0
150 real(kind_phys), intent(out ), dimension(:) :: sr
151
152 real(kind_phys), intent(in) :: dtp ! physics time step
153 logical, intent (in) :: hydrostatic, fast_mp_consv
154
155 logical, intent (in) :: lradar
156 real(kind=kind_phys), intent(inout), dimension(:,:) :: refl_10cm
157 logical, intent (in) :: reset, effr_in
158 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: rew, rei, rer, res, reg
159 logical, intent (in) :: cplchm
160 ! ice and liquid water 3d precipitation fluxes - only allocated if cplchm is .true.
161 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: pfi_lsan, pfl_lsan
162
163 character(len=*), intent(out) :: errmsg
164 integer, intent(out) :: errflg
165
166 ! local variables
167 integer :: iis, iie, jjs, jje, kks, kke, kbot, ktop
168 integer :: i, k, kk
169 real(kind=kind_phys), dimension(1:im,1:levs) :: delp, dz, uin, vin, pt, qv1, ql1, qi1, qr1, qs1, qg1, &
170 qa1, qnl, qni, pt_dt, qa_dt, u_dt, v_dt, w, qv_dt, ql_dt,&
171 qr_dt, qi_dt, qs_dt, qg_dt, p123, refl
172 real(kind=kind_phys), dimension(1:im,1:levs) :: q_con, cappa !for inline MP option
173 real(kind=kind_phys), dimension(1:im,1,1:levs) :: pfils, pflls
174 real(kind=kind_phys), dimension(1:im,1,1:levs) :: adj_vmr, te
175 real(kind=kind_phys), dimension(1:im,1:levs) :: prefluxw, prefluxr, prefluxi, prefluxs, prefluxg
176 real(kind=kind_phys), dimension(1:im) :: hs, gsize
177 real(kind=kind_dbl_prec), dimension(1:im) :: dte
178 !real(kind=kind_phys), dimension(:,:), allocatable :: den
179 real(kind=kind_phys), dimension(1:im) :: water0
180 real(kind=kind_phys) :: onebg
181 real(kind=kind_phys) :: tem
182 logical last_step, do_inline_mp
183
184 ! Initialize CCPP error handling variables
185 errmsg = ''
186 errflg = 0
187
188 iis = 1
189 iie = im
190 jjs = 1
191 jje = 1
192 kks = 1
193 kke = levs
194 ! flipping of vertical direction
195 ktop = 1
196 kbot = levs
197
198 onebg = con_one/con_g
199
200 do k = 1, levs
201 kk = levs-k+1
202 do i = 1, im
203 qv_dt(i,k) = 0.0
204 ql_dt(i,k) = 0.0
205 qr_dt(i,k) = 0.0
206 qi_dt(i,k) = 0.0
207 qs_dt(i,k) = 0.0
208 qg_dt(i,k) = 0.0
209 qa_dt(i,k) = 0.0
210 pt_dt(i,k) = 0.0
211 u_dt(i,k) = 0.0
212 v_dt(i,k) = 0.0
213 qnl(i,k) = aerfld(i,kk,11) ! sulfate
214 pfils(i,1,k) = 0.0
215 pflls(i,1,k) = 0.0
216 prefluxw(i,k) =0.0
217 prefluxi(i,k) =0.0
218 prefluxr(i,k) =0.0
219 prefluxs(i,k) =0.0
220 prefluxg(i,k) =0.0
221
222 ! flip vertical (k) coordinate top =1
223 qv1(i,k) = gq0(i,kk)
224 ql1(i,k) = gq0_ntcw(i,kk)
225 qr1(i,k) = gq0_ntrw(i,kk)
226 qi1(i,k) = gq0_ntiw(i,kk)
227 qs1(i,k) = gq0_ntsw(i,kk)
228 qg1(i,k) = gq0_ntgl(i,kk)
229 qa1(i,k) = gq0_ntclamt(i,kk)
230 pt(i,k) = gt0(i,kk)
231 w(i,k) = -vvl(i,kk) * (con_one+con_fvirt * gq0(i,kk)) &
232 * gt0(i,kk) / prsl(i,kk) * (con_rd*onebg)
233 uin(i,k) = gu0(i,kk)
234 vin(i,k) = gv0(i,kk)
235 delp(i,k) = del(i,kk)
236 dz(i,k) = (phii(i,kk)-phii(i,kk+1))*onebg
237 p123(i,k) = prsl(i,kk)
238 qni(i,k) = 10.
239 q_con(i,k) = 0.0
240 cappa(i,k) = 0.0
241 enddo
242 enddo
243
244 ! reset precipitation amounts to zero
245 water0 = 0
246 rain0 = 0
247 ice0 = 0
248 snow0 = 0
249 graupel0 = 0
250
251 ! Call MP driver
252 last_step = .false.
253 do_inline_mp = .false.
254 hs = oro(:) * con_g
255 gsize = sqrt(garea(:))
256
257 call gfdl_cloud_microphys_v3_mod_driver( qv1, ql1, qr1, qi1, qs1, qg1, qa1, qnl, qni, pt, w,&
258 uin, vin, dz, delp, gsize, dtp, hs, water0, rain0, &
259 ice0, snow0, graupel0, hydrostatic, iis, iie, kks, kke, q_con, cappa, &
260 fast_mp_consv, adj_vmr, te, dte, prefluxw, prefluxr, prefluxi, prefluxs, &
261 prefluxg, last_step, do_inline_mp )
262 tem = dtp*con_p001/con_secinday
263
264 ! fix negative values
265 do i = 1, im
266 !rain0(i) = max(con_d00, rain0(i))
267 !snow0(i) = max(con_d00, snow0(i))
268 !ice0(i) = max(con_d00, ice0(i))
269 !graupel0(i) = max(con_d00, graupel0(i))
270 if(water0(i)*tem < rainmin) then
271 water0(i) = 0.0
272 endif
273 if(rain0(i)*tem < rainmin) then
274 rain0(i) = 0.0
275 endif
276 if(ice0(i)*tem < rainmin) then
277 ice0(i) = 0.0
278 endif
279 if(snow0(i)*tem < rainmin) then
280 snow0(i) = 0.0
281 endif
282 if(graupel0(i)*tem < rainmin) then
283 graupel0(i) = 0.0
284 endif
285 enddo
286
287 ! calculate fraction of frozen precipitation using unscaled
288 ! values of rain0, ice0, snow0, graupel0 (for bit-for-bit)
289 do i=1,im
290 prcp0(i) = (rain0(i)+snow0(i)+ice0(i)+graupel0(i)) * tem
291 if ( prcp0(i) > rainmin ) then
292 sr(i) = (snow0(i) + ice0(i) + graupel0(i)) &
293 / (rain0(i) + snow0(i) + ice0(i) + graupel0(i))
294 else
295 sr(i) = 0.0
296 endif
297 enddo
298
299 ! convert rain0, ice0, snow0, graupel0 from mm per day to m per physics timestep
300 water0 = water0*tem
301 rain0 = rain0*tem
302 ice0 = ice0*tem
303 snow0 = snow0*tem
304 graupel0 = graupel0*tem
305
306 ! flip vertical coordinate back
307 do k=1,levs
308 kk = levs-k+1
309 do i=1,im
310 gq0(i,k) = qv1(i,kk)
311 gq0_ntcw(i,k) = ql1(i,kk)
312 gq0_ntrw(i,k) = qr1(i,kk)
313 gq0_ntiw(i,k) = qi1(i,kk)
314 gq0_ntsw(i,k) = qs1(i,kk)
315 gq0_ntgl(i,k) = qg1(i,kk)
316 gq0_ntclamt(i,k) = qa1(i,kk)
317 gt0(i,k) = pt(i,kk)
318 gu0(i,k) = uin(i,kk)
319 gv0(i,k) = vin(i,kk)
320 refl_10cm(i,k) = refl(i,kk)
321 enddo
322 enddo
323
324 ! output ice and liquid water 3d precipitation fluxes if requested
325 if (cplchm) then
326 do k=1,levs
327 kk = levs-k+1
328 do i=1,im
329 pfi_lsan(i,k) = prefluxi(i,kk) + prefluxs(i,kk) + prefluxg(i,kk)
330 pfl_lsan(i,k) = prefluxr(i,kk)
331 enddo
332 enddo
333 endif
334
335 if(effr_in) then
336 call cld_eff_rad (1, im, 1, levs, slmsk(1:im), &
337 prsl(1:im,1:levs), del(1:im,1:levs), &
338 gt0(1:im,1:levs), gq0(1:im,1:levs), &
339 gq0_ntcw(1:im,1:levs), gq0_ntiw(1:im,1:levs), &
340 gq0_ntrw(1:im,1:levs), gq0_ntsw(1:im,1:levs), &
341 gq0_ntgl(1:im,1:levs), gq0_ntclamt(1:im,1:levs), &
342 rew(1:im,1:levs), rei(1:im,1:levs), rer(1:im,1:levs),&
343 res(1:im,1:levs), reg(1:im,1:levs),snowd(1:im))
344 endif
345
346 if(lradar) then
347 call rad_ref (1, im, 1, 1, qv1(1:im,1:levs), qr1(1:im,1:levs), &
348 qs1(1:im,1:levs),qg1(1:im,1:levs),pt(1:im,1:levs), &
349 delp(1:im,1:levs), dz(1:im,1:levs), refl(1:im,1:levs), levs, hydrostatic, &
350 do_inline_mp, 1)
351
352 do k=1,levs
353 kk = levs-k+1
354 do i=1,im
355 refl_10cm(i,k) = max(-35.,refl(i,kk))
356 enddo
357 enddo
358 endif
359
360 end subroutine gfdl_cloud_microphys_v3_run
361