CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
gfdl_cloud_microphys.F90
1
4
8
13
14 implicit none
15
16 private
17
18 public gfdl_cloud_microphys_run, gfdl_cloud_microphys_init, gfdl_cloud_microphys_finalize
19
20 logical :: is_initialized = .false.
21
22contains
23
24! -----------------------------------------------------------------------
25! CCPP entry points for gfdl cloud microphysics
26! -----------------------------------------------------------------------
27
34 subroutine gfdl_cloud_microphys_init (me, master, nlunit, input_nml_file, logunit, fn_nml, &
35 imp_physics, imp_physics_gfdl, do_shoc, errmsg, errflg)
36
37 implicit none
38
39 integer, intent (in) :: me
40 integer, intent (in) :: master
41 integer, intent (in) :: nlunit
42 integer, intent (in) :: logunit
43 character(len=*), intent (in) :: fn_nml
44 character(len=*), intent (in) :: input_nml_file(:)
45 integer, intent( in) :: imp_physics
46 integer, intent( in) :: imp_physics_gfdl
47 logical, intent( in) :: do_shoc
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'
65 errflg = 1
66 return
67 endif
68
69 call gfdl_cloud_microphys_mod_init(me, master, nlunit, input_nml_file, logunit, fn_nml, errmsg, errflg)
70
71 is_initialized = .true.
72
73 end subroutine gfdl_cloud_microphys_init
74
75! =======================================================================
82 subroutine gfdl_cloud_microphys_finalize(errmsg, errflg)
83
84 implicit none
85
86 character(len=*), intent(out) :: errmsg
87 integer, intent(out) :: errflg
88
89 ! Initialize CCPP error handling variables
90 errmsg = ''
91 errflg = 0
92
93 if (.not.is_initialized) return
94
96
97 is_initialized = .false.
98
99 end subroutine gfdl_cloud_microphys_finalize
100
117 subroutine gfdl_cloud_microphys_run( &
118 levs, im, rainmin, con_g, con_fvirt, con_rd, con_eps, frland, garea, islmsk, &
119 gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, gq0_ntsw, gq0_ntgl, gq0_ntclamt, &
120 gt0, gu0, gv0, vvl, prsl, phii, del, &
121 rain0, ice0, snow0, graupel0, prcp0, sr, &
122 dtp, hydrostatic, phys_hydrostatic, lradar, refl_10cm, &
123 reset, effr_in, rew, rei, rer, res, reg, &
124 cplchm, pfi_lsan, pfl_lsan, errmsg, errflg)
125
126 use machine, only: kind_phys
127
128 implicit none
129
130 ! DH* TODO: CLEANUP, all of these should be coming in through the argument list
131 ! parameters
132 real(kind=kind_phys), parameter :: one = 1.0d0
133 real(kind=kind_phys), parameter :: con_p001= 0.001d0
134 real(kind=kind_phys), parameter :: con_day = 86400.d0
135 !real(kind=kind_phys), parameter :: rainmin = 1.0d-13
136 ! *DH
137
138 ! interface variables
139 integer, intent(in ) :: levs, im
140 real(kind=kind_phys), intent(in ) :: con_g, con_fvirt, con_rd, con_eps, rainmin
141 real(kind=kind_phys), intent(in ), dimension(:) :: frland, garea
142 integer, intent(in ), dimension(:) :: islmsk
143 real(kind=kind_phys), intent(inout), dimension(:,:) :: gq0, gq0_ntcw, gq0_ntrw, gq0_ntiw, &
144 gq0_ntsw, gq0_ntgl, gq0_ntclamt
145 real(kind=kind_phys), intent(inout), dimension(:,:) :: gt0, gu0, gv0
146 real(kind=kind_phys), intent(in ), dimension(:,:) :: vvl, prsl, del
147 real(kind=kind_phys), intent(in ), dimension(:,:) :: phii
148
149 ! rain/snow/ice/graupel/precip amounts, fraction of frozen precip
150 real(kind_phys), intent(out ), dimension(:), optional :: rain0
151 real(kind_phys), intent(out ), dimension(:), optional :: snow0
152 real(kind_phys), intent(out ), dimension(:), optional :: ice0
153 real(kind_phys), intent(out ), dimension(:), optional :: graupel0
154 real(kind_phys), intent(out ), dimension(:) :: prcp0
155 real(kind_phys), intent(out ), dimension(:) :: sr
156
157 real(kind_phys), intent(in) :: dtp ! physics time step
158 logical, intent (in) :: hydrostatic, phys_hydrostatic
159
160 logical, intent (in) :: lradar
161 real(kind=kind_phys), intent(inout), dimension(:,:) :: refl_10cm
162 logical, intent (in) :: reset, effr_in
163 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: rew, rei, rer, res, reg
164 logical, intent (in) :: cplchm
165 ! ice and liquid water 3d precipitation fluxes - only allocated if cplchm is .true.
166 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: pfi_lsan, pfl_lsan
167
168 character(len=*), intent(out) :: errmsg
169 integer, intent(out) :: errflg
170
171 ! local variables
172 integer :: iis, iie, jjs, jje, kks, kke, kbot, ktop
173 integer :: i, k, kk
174 real(kind=kind_phys), dimension(1:im,1:levs) :: delp, dz, uin, vin, pt, qv1, ql1, qr1, qg1, qa1, qn1, qi1, &
175 qs1, pt_dt, qa_dt, u_dt, v_dt, w, qv_dt, ql_dt, qr_dt, qi_dt, &
176 qs_dt, qg_dt, p123, refl
177 real(kind=kind_phys), dimension(1:im,1,1:levs) :: pfils, pflls
178 real(kind=kind_phys), dimension(:,:), allocatable :: den
179 real(kind=kind_phys) :: onebg
180 real(kind=kind_phys) :: tem
181
182 ! Initialize CCPP error handling variables
183 errmsg = ''
184 errflg = 0
185
186 iis = 1
187 iie = im
188 jjs = 1
189 jje = 1
190 kks = 1
191 kke = levs
192 ! flipping of vertical direction
193 ktop = 1
194 kbot = levs
195
196 onebg = one/con_g
197
198 do k = 1, levs
199 kk = levs-k+1
200 do i = 1, im
201 qv_dt(i,k) = 0.0
202 ql_dt(i,k) = 0.0
203 qr_dt(i,k) = 0.0
204 qi_dt(i,k) = 0.0
205 qs_dt(i,k) = 0.0
206 qg_dt(i,k) = 0.0
207 qa_dt(i,k) = 0.0
208 pt_dt(i,k) = 0.0
209 u_dt(i,k) = 0.0
210 v_dt(i,k) = 0.0
211 qn1(i,k) = 0.0
212 pfils(i,1,k) = 0.0
213 pflls(i,1,k) = 0.0
214 ! flip vertical (k) coordinate
215 qv1(i,k) = gq0(i,kk)
216 ql1(i,k) = gq0_ntcw(i,kk)
217 qr1(i,k) = gq0_ntrw(i,kk)
218 qi1(i,k) = gq0_ntiw(i,kk)
219 qs1(i,k) = gq0_ntsw(i,kk)
220 qg1(i,k) = gq0_ntgl(i,kk)
221 qa1(i,k) = gq0_ntclamt(i,kk)
222 pt(i,k) = gt0(i,kk)
223 w(i,k) = -vvl(i,kk) * (one+con_fvirt * gq0(i,kk)) &
224 * gt0(i,kk) / prsl(i,kk) * (con_rd*onebg)
225 uin(i,k) = gu0(i,kk)
226 vin(i,k) = gv0(i,kk)
227 delp(i,k) = del(i,kk)
228 dz(i,k) = (phii(i,kk)-phii(i,kk+1))*onebg
229 p123(i,k) = prsl(i,kk)
230 refl(i,k) = refl_10cm(i,kk)
231 enddo
232 enddo
233
234 ! reset precipitation amounts to zero
235 rain0 = 0
236 ice0 = 0
237 snow0 = 0
238 graupel0 = 0
239
240 call gfdl_cloud_microphys_mod_driver(iis, iie, jjs, jje, kks, kke, ktop, kbot, &
241 qv1, ql1, qr1, qi1, qs1, qg1, qa1, qn1, qv_dt, ql_dt, qr_dt, qi_dt, &
242 qs_dt, qg_dt, qa_dt, pt_dt, pt, w, uin, vin, u_dt, v_dt, dz, delp, &
243 garea, dtp, frland, rain0, snow0, ice0, graupel0, hydrostatic, &
244 phys_hydrostatic, p123, lradar, refl, reset, pfils, pflls)
245 tem = dtp*con_p001/con_day
246
247 ! fix negative values
248 do i = 1, im
249 !rain0(i) = max(con_d00, rain0(i))
250 !snow0(i) = max(con_d00, snow0(i))
251 !ice0(i) = max(con_d00, ice0(i))
252 !graupel0(i) = max(con_d00, graupel0(i))
253 if(rain0(i)*tem < rainmin) then
254 rain0(i) = 0.0
255 endif
256 if(ice0(i)*tem < rainmin) then
257 ice0(i) = 0.0
258 endif
259 if(snow0(i)*tem < rainmin) then
260 snow0(i) = 0.0
261 endif
262 if(graupel0(i)*tem < rainmin) then
263 graupel0(i) = 0.0
264 endif
265 enddo
266
267 ! calculate fraction of frozen precipitation using unscaled
268 ! values of rain0, ice0, snow0, graupel0 (for bit-for-bit)
269 do i=1,im
270 prcp0(i) = (rain0(i)+snow0(i)+ice0(i)+graupel0(i)) * tem
271 if ( prcp0(i) > rainmin ) then
272 sr(i) = (snow0(i) + ice0(i) + graupel0(i)) &
273 / (rain0(i) + snow0(i) + ice0(i) + graupel0(i))
274 else
275 sr(i) = 0.0
276 endif
277 enddo
278
279 ! convert rain0, ice0, snow0, graupel0 from mm per day to m per physics timestep
280 rain0 = rain0*tem
281 ice0 = ice0*tem
282 snow0 = snow0*tem
283 graupel0 = graupel0*tem
284
285 ! flip vertical coordinate back
286 do k=1,levs
287 kk = levs-k+1
288 do i=1,im
289 gq0(i,k) = qv1(i,kk) + qv_dt(i,kk) * dtp
290 gq0_ntcw(i,k) = ql1(i,kk) + ql_dt(i,kk) * dtp
291 gq0_ntrw(i,k) = qr1(i,kk) + qr_dt(i,kk) * dtp
292 gq0_ntiw(i,k) = qi1(i,kk) + qi_dt(i,kk) * dtp
293 gq0_ntsw(i,k) = qs1(i,kk) + qs_dt(i,kk) * dtp
294 gq0_ntgl(i,k) = qg1(i,kk) + qg_dt(i,kk) * dtp
295 gq0_ntclamt(i,k) = qa1(i,kk) + qa_dt(i,kk) * dtp
296 gt0(i,k) = gt0(i,k) + pt_dt(i,kk) * dtp
297 gu0(i,k) = gu0(i,k) + u_dt(i,kk) * dtp
298 gv0(i,k) = gv0(i,k) + v_dt(i,kk) * dtp
299 refl_10cm(i,k) = refl(i,kk)
300 enddo
301 enddo
302
303 ! output ice and liquid water 3d precipitation fluxes if requested
304 if (cplchm) then
305 do k=1,levs
306 kk = levs-k+1
307 do i=1,im
308 pfi_lsan(i,k) = pfils(i,1,kk)
309 pfl_lsan(i,k) = pflls(i,1,kk)
310 enddo
311 enddo
312 endif
313
314 if(effr_in) then
315 allocate(den(1:im,1:levs))
316 do k=1,levs
317 do i=1,im
318 den(i,k)=con_eps*prsl(i,k)/(con_rd*gt0(i,k)*(gq0(i,k)+con_eps))
319 enddo
320 enddo
321 call cloud_diagnosis (1, im, 1, levs, den(1:im,1:levs), &
322 del(1:im,1:levs), islmsk(1:im), &
323 gq0_ntcw(1:im,1:levs), gq0_ntiw(1:im,1:levs), &
324 gq0_ntrw(1:im,1:levs), &
325 gq0_ntsw(1:im,1:levs) + gq0_ntgl(1:im,1:levs), &
326 gq0_ntgl(1:im,1:levs)*0.0, gt0(1:im,1:levs), &
327 rew(1:im,1:levs), rei(1:im,1:levs), rer(1:im,1:levs),&
328 res(1:im,1:levs), reg(1:im,1:levs))
329 deallocate(den)
330 endif
331
332 end subroutine gfdl_cloud_microphys_run
333
334end module gfdl_cloud_microphys
subroutine, public cloud_diagnosis(is, ie, ks, ke, den, delp, lsm, qmw, qmi, qmr, qms, qmg, t, rew, rei, rer, res, reg)
The subroutine 'cloud_diagnosis' diagnoses the radius of cloud species.
subroutine, public gfdl_cloud_microphys_mod_end()
The subroutine 'gfdl_cloud_microphys_init' terminates the GFDL cloud microphysics.
subroutine, public gfdl_cloud_microphys_mod_driver(iis, iie, jjs, jje, kks, kke, ktop, kbot, qv, ql, qr, qi, qs, qg, qa, qn, qv_dt, ql_dt, qr_dt, qi_dt, qs_dt, qg_dt, qa_dt, pt_dt, pt, w, uin, vin, udt, vdt, dz, delp, area, dt_in, land, rain, snow, ice, graupel, hydrostatic, phys_hydrostatic, p, lradar, refl_10cm, reset, pfils, pflls)
This subroutine is the driver of the GFDL cloud microphysics.
subroutine, public gfdl_cloud_microphys_mod_init(me, master, nlunit, input_nml_file, logunit, fn_nml, errmsg, errflg)
The subroutine 'gfdl_cloud_microphys_init' initializes the GFDL cloud microphysics.
This module contains the column GFDL Cloud microphysics scheme.
This module contains the CCPP entry point for the column GFDL cloud microphysics ( Chen and Lin (2013...