Line data Source code
1 : ! Compute tendencies from sedimentation of cloud liquid and ice particles
2 : ! Original authors: Byron Boville, September 2002 from code by P.J. Rasch
3 : ! CCPP-ized: Haipeng Lin, January 2025
4 : module cloud_particle_sedimentation
5 :
6 : use ccpp_kinds, only: kind_phys
7 :
8 : implicit none
9 : private
10 : save
11 :
12 : ! public CCPP-compliant subroutines
13 : public :: cloud_particle_sedimentation_init
14 : public :: cloud_particle_sedimentation_run
15 :
16 : ! tuning parameters for cloud liquid and ice particle sedimentation
17 : real(kind_phys), parameter :: vland = 1.5_kind_phys ! liquid fall velocity over land [cm s-1]
18 : real(kind_phys), parameter :: vocean = 2.8_kind_phys ! liquid fall velocity over ocean [cm s-1]
19 : real(kind_phys), parameter :: mxsedfac = 0.99_kind_phys ! maximum sedimentation flux factor
20 :
21 : ! use Stokes velocity method
22 : ! or McFarquhar and Heymsfield (https://doi.org/10.1175/1520-0469(1997)054<2187:POTCIC>2.0.CO;2)
23 : logical, parameter :: stokes = .true. ! use Stokes velocity instead of McFarquhar and Heymsfield
24 :
25 : ! tuning parameters for Stokes terminal velocity
26 : real(kind_phys) :: cldsed_ice_stokes_fac ! factor applied to ice fall vel
27 : ! from Stokes terminal velocity
28 : real(kind_phys), parameter :: eta = 1.7e-5_kind_phys ! viscosity of air [kg m s-1]
29 : real(kind_phys), parameter :: r40 = 40._kind_phys ! 40 micron radius
30 : real(kind_phys), parameter :: r400 = 400._kind_phys ! 400 micron radius
31 : real(kind_phys), parameter :: v400 = 1.00_kind_phys ! fall velocity of 400 micron sphere [m s-1]
32 : real(kind_phys) :: v40 ! Stokes fall velocity of 40 micron sphere [m s-1]
33 : real(kind_phys) :: vslope ! linear slope for large particles [m s-1 micron-1]
34 :
35 : ! tuning parameters for McFarquhar and Heymsfield terminal velocity
36 : real(kind_phys), parameter :: vice_small = 1._kind_phys ! ice fall velocity for small concentration [cm s-1]
37 : real(kind_phys) :: lbound ! lower bound for iciwc
38 :
39 : ! misc
40 : real(kind_phys), parameter :: um_to_m = 1.e-12_kind_phys
41 :
42 : contains
43 :
44 : !> \section arg_table_cloud_particle_sedimentation_init Argument Table
45 : !! \htmlinclude arg_table_cloud_particle_sedimentation_init.html
46 1024 : subroutine cloud_particle_sedimentation_init(&
47 : amIRoot, iulog, &
48 : cldsed_ice_stokes_fac_in, &
49 : rhoh2o, gravit, &
50 1024 : errmsg, errflg)
51 :
52 : ! Input arguments
53 : logical, intent(in) :: amIRoot
54 : integer, intent(in) :: iulog ! log output unit
55 : real(kind_phys), intent(in) :: cldsed_ice_stokes_fac_in
56 : real(kind_phys), intent(in) :: gravit
57 : real(kind_phys), intent(in) :: rhoh2o
58 :
59 : ! Output arguments
60 : character(len=512), intent(out) :: errmsg ! error message
61 : integer, intent(out) :: errflg ! error flag
62 :
63 : ! Local variables
64 : real(kind_phys) :: ac, bc, cc ! constants for McFarquhar and Heymsfield method
65 :
66 1024 : errmsg = ''
67 1024 : errflg = 0
68 1024 : cldsed_ice_stokes_fac = cldsed_ice_stokes_fac_in
69 :
70 : ! linear ramp variables for fall velocity
71 1024 : v40 = (2._kind_phys/9._kind_phys)*rhoh2o*gravit/eta*r40**2*1.e-12_kind_phys
72 1024 : vslope = (v400 - v40)/(r400 - r40)
73 :
74 : ! McFarquhar and Heymsfield lower bound for iciwc
75 1024 : cc = 128.64_kind_phys
76 1024 : bc = 53.242_kind_phys
77 1024 : ac = 5.4795_kind_phys
78 1024 : lbound = (-bc + sqrt(bc*bc - 4*ac*cc))/(2*ac)
79 1024 : lbound = 10._kind_phys**lbound
80 :
81 1024 : if(amIRoot) then
82 2 : write(iulog,*) 'cloud_particle_sedimentation_init: cldsed_ice_stokes_fac = ', cldsed_ice_stokes_fac
83 : endif
84 :
85 1024 : end subroutine cloud_particle_sedimentation_init
86 :
87 : ! Compute gravitational sedimentation velocities for cloud liquid water
88 : ! and ice, based on Lawrence and Crutzen (1998).
89 : ! https://doi.org/10.3402/tellusb.v50i3.16129
90 : !
91 : ! and apply cloud particle gravitational sedimentation to condensate.
92 : !
93 : ! Original authors: mgl, March 1998 (MATCH-MPIC version 2.0)
94 : ! adapted by P.J. Rasch
95 : ! B.A. Boville, September 2002
96 : ! P.J. Rasch, May 2003
97 : !> \section arg_table_cloud_particle_sedimentation_run Argument Table
98 : !! \htmlinclude arg_table_cloud_particle_sedimentation_run.html
99 70392 : subroutine cloud_particle_sedimentation_run( &
100 : ncol, pver, pverp, dtime, &
101 : tmelt, gravit, latvap, latice, rair, rhoh2o, &
102 : icritc, & ! from prognostic_cloud_water namelist
103 211176 : pint, pmid, pdel, t, cloud, &
104 140784 : icefrac, landfrac, ocnfrac, &
105 281568 : cldliq, cldice, snowh, landm, &
106 : ! below output for sedimentation velocity
107 140784 : pvliq, pvice, &
108 : ! below output for sedimentation to condensate
109 70392 : liqtend, icetend, wvtend, htend, sfliq, sfice, &
110 0 : errmsg, errflg)
111 : ! To-be-CCPPized dependencies
112 : use cloud_optical_properties, only: reltab, reitab
113 :
114 : ! Input arguments
115 : integer, intent(in) :: ncol
116 : integer, intent(in) :: pver
117 : integer, intent(in) :: pverp
118 : real(kind_phys), intent(in) :: dtime
119 : real(kind_phys), intent(in) :: tmelt
120 : real(kind_phys), intent(in) :: gravit
121 : real(kind_phys), intent(in) :: latvap
122 : real(kind_phys), intent(in) :: latice
123 : real(kind_phys), intent(in) :: rair
124 : real(kind_phys), intent(in) :: rhoh2o
125 : real(kind_phys), intent(in) :: icritc ! tunable_parameter_for_autoconversion_of_cold_ice_for_rk_microphysics [kg kg-1]
126 : real(kind_phys), intent(in) :: pint(:,:) ! air_pressure_at_interface [Pa]
127 : real(kind_phys), intent(in) :: pmid(:,:) ! air_pressure [Pa]
128 : real(kind_phys), intent(in) :: pdel(:,:) ! air_pressure_thickness [Pa]
129 : real(kind_phys), intent(in) :: t(:,:) ! air_temperature [K]
130 : real(kind_phys), intent(in) :: cloud(:,:) ! cloud_area_fraction [fraction]
131 : real(kind_phys), intent(in) :: icefrac(:) ! sea_ice_area_fraction [fraction]
132 : real(kind_phys), intent(in) :: landfrac(:) ! land_area_fraction [fraction]
133 : real(kind_phys), intent(in) :: ocnfrac(:) ! ocean_area_fraction [fraction]
134 : real(kind_phys), intent(in) :: cldliq(:,:) ! adv: cloud_liquid_water_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
135 : real(kind_phys), intent(in) :: cldice(:,:) ! adv: cloud_ice_mixing_ratio_wrt_moist_air_and_condensed_water [kg kg-1]
136 : real(kind_phys), intent(in) :: snowh(:) ! lwe_surface_snow_depth_over_land [m]
137 : real(kind_phys), intent(in) :: landm(:) ! smoothed_land_area_fraction [fraction]
138 :
139 : ! Output arguments
140 : ! note: pvliq, pvice are at the interfaces (loss from cell is based on pvel(k+1))
141 : real(kind_phys), intent(out) :: pvliq(:,:) ! vertical velocity of cloud liquid drops [Pa s-1]
142 : real(kind_phys), intent(out) :: pvice(:,:) ! vertical velocity of cloud ice particles [Pa s-1]
143 :
144 : real(kind_phys), intent(out) :: liqtend(:,:) ! liquid condensate tendency [kg kg-1 s-1] -- to apply cldliq tendency
145 : real(kind_phys), intent(out) :: icetend(:,:) ! ice condensate tendency [kg kg-1 s-1] -- to apply cldice tendency
146 : real(kind_phys), intent(out) :: wvtend(:,:) ! water vapor tendency [kg kg-1 s-1] -- to apply wv tendency
147 : real(kind_phys), intent(out) :: htend(:,:) ! heating rate [J kg-1 s-1] -- to apply s tendency
148 : real(kind_phys), intent(out) :: sfliq(:) ! surface flux of liquid (rain) [kg m-2 s-1]
149 : real(kind_phys), intent(out) :: sfice(:) ! lwe_cloud_ice_sedimentation_rate_at_surface_due_to_microphysics [m s-1]
150 : character(len=512), intent(out) :: errmsg ! error message
151 : integer, intent(out) :: errflg ! error flag
152 :
153 : ! Local variables
154 140784 : real(kind_phys) :: rho(ncol, pver) ! air density [kg m-3]
155 : real(kind_phys) :: vfall ! settling velocity of cloud particles [m s-1]
156 : real(kind_phys) :: icice ! in cloud ice water content [kg kg-1]
157 : real(kind_phys) :: iciwc ! in cloud ice water content [g m-3]
158 : real(kind_phys) :: icefac
159 : real(kind_phys) :: logiwc
160 140784 : real(kind_phys) :: rei(ncol, pver) ! effective radius of ice particles [um]
161 140784 : real(kind_phys) :: rel(ncol, pver) ! effective radius of liquid particles [um]
162 140784 : real(kind_phys) :: fxliq(ncol, pverp) ! fluxes at interface, liquid (positive is down) [Pa]
163 140784 : real(kind_phys) :: fxice(ncol, pverp) ! fluxes at interface, ice (positive is down) [Pa]
164 140784 : real(kind_phys) :: cldab(ncol) ! cloud in layer above [fraction]
165 : real(kind_phys) :: evapliq ! evaporation of cloud liquid into environment [kg kg-1 s-1]
166 : real(kind_phys) :: evapice ! evaporation of cloud ice into environment [kg kg-1 s-1]
167 : real(kind_phys) :: cldovrl ! cloud overlap factor
168 : integer :: i, k
169 :
170 70392 : errmsg = ''
171 70392 : errflg = 0
172 :
173 : !--------------------------------------------------
174 : ! COMPUTE GRAVITATIONAL SEDIMENTATION VELOCITIES
175 : !--------------------------------------------------
176 : ! NEED TO BE CAREFUL - VELOCITIES SHOULD BE AT THE *LOWER* INTERFACE
177 : ! (THAT IS, FOR K+1), FLUXES ARE ALSO AT THE LOWER INTERFACE (K+1),
178 : ! BUT MIXING RATIOS ARE AT THE MIDPOINTS (K)...
179 : !
180 : ! NOTE THAT PVEL IS ON INTERFACES, WHEREAS VFALL IS FOR THE CELL
181 : ! AVERAGES (I.E., MIDPOINTS); ASSUME THE FALL VELOCITY APPLICABLE TO THE
182 : ! LOWER INTERFACE (K+1) IS THE SAME AS THAT APPLICABLE FOR THE CELL (V(K))
183 :
184 : ! LIQUID:
185 : ! The fall velocities assume that droplets have a gamma distribution
186 : ! with effective radii for land and ocean as assessed by Han et al.;
187 : ! see Lawrence and Crutzen (1998) for a derivation.
188 28436184 : rho(:ncol,:) = pmid(:ncol,:)/(rair*t(:ncol,:))
189 29527176 : pvliq(:ncol,:) = 0._kind_phys
190 :
191 : ! get effective radius of liquid drop (rel)
192 70392 : call reltab(ncol=ncol, pver=pver, tmelt=tmelt, t=t(:ncol,:), landfrac=landfrac(:ncol), landm=landm(:ncol), icefrac=icefrac(:ncol), snowh=snowh(:ncol), rel=rel(:ncol,:))
193 :
194 1900584 : do k = 1, pver
195 28436184 : do i = 1, ncol
196 28365792 : if (cloud(i, k) > 0._kind_phys .and. cldliq(i, k) > 0._kind_phys) then
197 5651785 : if (rel(i, k) < 40._kind_phys) then
198 5651785 : vfall = 2._kind_phys/9._kind_phys*rhoh2o*gravit*rel(i, k)**2/eta*um_to_m ! um^2 -> m^2
199 : else
200 0 : vfall = v40 + vslope*(rel(i, k) - r40) ! linear above 40 microns
201 : end if
202 : ! convert the fall speed to pressure units, but do not apply the traditional
203 : ! negative convention for pvel.
204 5651785 : pvliq(i, k + 1) = vfall*rho(i, k)*gravit ! m s-1 to Pa s-1
205 : end if
206 : end do
207 : end do
208 :
209 : ! ICE:
210 : ! The fall velocities are based on data from McFarquhar and Heymsfield
211 : ! or on Stokes terminal velocity for spheres and the effective radius.
212 29527176 : pvice(:ncol,:) = 0._kind_phys
213 :
214 : if(stokes) then
215 : ! stokes terminal velocity < 40 microns
216 :
217 : ! get effective radius (rei)
218 70392 : call reitab(ncol=ncol, pver=pver, t=t(:ncol,:), re=rei(:ncol,:))
219 :
220 1900584 : do k = 1, pver
221 28436184 : do i = 1, ncol
222 28365792 : if (cloud(i, k) > 0._kind_phys .and. cldice(i, k) > 0._kind_phys) then
223 4823028 : if (rei(i, k) < 40._kind_phys) then
224 1344565 : vfall = 2._kind_phys/9._kind_phys*rhoh2o*gravit*rei(i, k)**2/eta*um_to_m ! microns^2 -> m^2
225 1344565 : vfall = vfall*cldsed_ice_stokes_fac
226 : else
227 3478463 : vfall = v40 + vslope*(rei(i, k) - r40) ! linear above 40 microns
228 : end if
229 :
230 : ! convert the fall speed to pressure units, but do not apply the traditional
231 : ! negative convention for pvel.
232 4823028 : pvice(i, k + 1) = vfall*rho(i, k)*gravit ! m s-1 to Pa s-1
233 : end if
234 : end do
235 : end do
236 :
237 : else
238 : ! McFarquhar and Heymsfield > icritc
239 : do k = 1, pver
240 : do i = 1, ncol
241 : if (cloud(i, k) > 0._kind_phys .and. cldice(i, k) > 0._kind_phys) then
242 :
243 : ! compute the in-cloud ice concentration (kg/kg)
244 : icice = cldice(i, k)/cloud(i, k)
245 :
246 : ! compute the ice water content in g/m3
247 : iciwc = icice*rho(i, k)*1.e3_kind_phys
248 :
249 : ! set the fall velocity (cm/s) to depend on the ice water content in g/m3,
250 : if (iciwc > lbound) then ! need this because of log10
251 : logiwc = log10(iciwc)
252 : ! Median -
253 : vfall = 128.64_kind_phys + 53.242_kind_phys*logiwc + 5.4795_kind_phys*logiwc**2
254 : ! Average -
255 : ! vfall = 122.63 + 44.111*logiwc + 4.2144*logiwc**2
256 : else
257 : vfall = 0._kind_phys
258 : end if
259 :
260 : ! set ice velocity to 1 cm/s if ice mixing ratio < icritc, ramp to value
261 : ! calculated above at 2*icritc
262 : if (icice <= icritc) then
263 : vfall = vice_small
264 : else if (icice < 2*icritc) then
265 : icefac = (icice - icritc)/icritc
266 : vfall = vice_small*(1._kind_phys - icefac) + vfall*icefac
267 : end if
268 :
269 : ! bound the terminal velocity of ice particles at high concentration
270 : vfall = min(100.0_kind_phys, vfall)
271 :
272 : ! convert the fall speed to pressure units, but do not apply the traditional
273 : ! negative convention for pvel.
274 : pvice(i, k + 1) = vfall &
275 : *0.01_kind_phys & ! cm to meters
276 : *rho(i, k)*gravit ! m s-1 to Pa s-1
277 : end if
278 : end do
279 : end do
280 : endif
281 :
282 : !--------------------------------------------------
283 : ! APPLY CLOUD PARTICLE GRAVITATIONAL SEDIMENTATION TO CONDENSATE
284 : !--------------------------------------------------
285 : ! Initialize variables
286 29527176 : fxliq(:ncol, :) = 0._kind_phys ! flux at interfaces (liquid)
287 29527176 : fxice(:ncol, :) = 0._kind_phys ! flux at interfaces (ice)
288 28436184 : liqtend(:ncol, :) = 0._kind_phys ! condensate tend (liquid)
289 28436184 : icetend(:ncol, :) = 0._kind_phys ! condensate tend (ice)
290 28436184 : wvtend(:ncol, :) = 0._kind_phys ! environmental moistening
291 28436184 : htend(:ncol, :) = 0._kind_phys ! evaporative cooling
292 1090992 : sfliq(:ncol) = 0._kind_phys ! condensate sedimentation flux out bot of column (liquid)
293 1090992 : sfice(:ncol) = 0._kind_phys ! condensate sedimentation flux out bot of column (ice)
294 :
295 : ! Get fluxes at interior points
296 70392 : call getflx(ncol, pver, pverp, pint, cldliq, pvliq, dtime, fxliq, errmsg, errflg)
297 70392 : if(errflg /= 0) then
298 0 : return
299 : endif
300 70392 : call getflx(ncol, pver, pverp, pint, cldice, pvice, dtime, fxice, errmsg, errflg)
301 70392 : if(errflg /= 0) then
302 0 : return
303 : endif
304 :
305 : ! Calculate fluxes at boundaries
306 1090992 : do i = 1, ncol
307 1020600 : fxliq(i, 1) = 0._kind_phys
308 1020600 : fxice(i, 1) = 0._kind_phys
309 : ! surface flux by upstream scheme
310 1020600 : fxliq(i, pverp) = cldliq(i, pver)*pvliq(i, pverp)*dtime
311 1090992 : fxice(i, pverp) = cldice(i, pver)*pvice(i, pverp)*dtime
312 : end do
313 :
314 : ! filter out any negative fluxes from the getflx routine
315 : ! (typical fluxes are of order > 1e-3 when clouds are present)
316 1830192 : do k = 2, pver
317 27274800 : fxliq(:ncol, k) = max(0._kind_phys, fxliq(:ncol, k))
318 27345192 : fxice(:ncol, k) = max(0._kind_phys, fxice(:ncol, k))
319 : end do
320 :
321 : ! Limit the flux out of the bottom of each cell to the water content in each phase.
322 : ! Apply mxsedfac to prevent generating very small negative cloud water/ice
323 : ! NOTE, REMOVED CLOUD FACTOR FROM AVAILABLE WATER. ALL CLOUD WATER IS IN CLOUDS.
324 : ! ***Should we include the flux in the top, to allow for thin surface layers?
325 : ! ***Requires simple treatment of cloud overlap, already included below.
326 1900584 : do k = 1, pver
327 28436184 : do i = 1, ncol
328 26535600 : fxliq(i, k + 1) = min(fxliq(i, k + 1), mxsedfac*cldliq(i, k)*pdel(i, k))
329 28365792 : fxice(i, k + 1) = min(fxice(i, k + 1), mxsedfac*cldice(i, k)*pdel(i, k))
330 : ! fxliq(i,k+1) = min( fxliq(i,k+1), cldliq(i,k) * pdel(i,k) + fxliq(i,k))
331 : ! fxice(i,k+1) = min( fxice(i,k+1), cldice(i,k) * pdel(i,k) + fxice(i,k))
332 : ! fxliq(i,k+1) = min( fxliq(i,k+1), cloud(i,k) * cldliq(i,k) * pdel(i,k) )
333 : ! fxice(i,k+1) = min( fxice(i,k+1), cloud(i,k) * cldice(i,k) * pdel(i,k) )
334 : end do
335 : end do
336 :
337 : ! Now calculate the tendencies assuming that condensate evaporates when
338 : ! it falls into environment, and does not when it falls into cloud.
339 : ! All flux out of the layer comes from the cloudy part.
340 : ! Assume maximum overlap for stratiform clouds
341 : ! if cloud above < cloud, all water falls into cloud below
342 : ! if cloud above > cloud, water split between cloud and environment
343 1900584 : do k = 1, pver
344 28365792 : cldab(:ncol) = 0._kind_phys
345 28436184 : do i = 1, ncol
346 : ! cloud overlap cloud factor
347 26535600 : cldovrl = min(cloud(i, k)/(cldab(i) + .0001_kind_phys), 1._kind_phys)
348 26535600 : cldab(i) = cloud(i, k)
349 : ! evaporation into environment cause moistening and cooling
350 26535600 : evapliq = fxliq(i, k)*(1._kind_phys - cldovrl)/(dtime*pdel(i, k)) ! into env [kg kg-1 s-1]
351 26535600 : evapice = fxice(i, k)*(1._kind_phys - cldovrl)/(dtime*pdel(i, k)) ! into env [kg kg-1 s-1]
352 26535600 : wvtend(i, k) = evapliq + evapice ! evaporation into environment [kg kg-1 s-1]
353 26535600 : htend(i, k) = -latvap*evapliq - (latvap + latice)*evapice ! evaporation [W kg-1]
354 : ! net flux into cloud changes cloud liquid/ice (all flux is out of cloud)
355 26535600 : liqtend(i, k) = (fxliq(i, k)*cldovrl - fxliq(i, k + 1))/(dtime*pdel(i, k))
356 28365792 : icetend(i, k) = (fxice(i, k)*cldovrl - fxice(i, k + 1))/(dtime*pdel(i, k))
357 : end do
358 : end do
359 :
360 : ! convert flux out the bottom to mass units Pa -> kg/m2/s
361 1090992 : sfliq(:ncol) = fxliq(:ncol, pverp)/(dtime*gravit)
362 1090992 : sfice(:ncol) = fxice(:ncol, pverp)/(dtime*gravit)
363 :
364 : ! Convert lwe_cloud_ice_sedimentation_rate_at_surface_due_to_microphysics from kg m-2 s-1 to precip units m s-1
365 1090992 : sfice(:ncol) = sfice(:ncol)/1000._kind_phys
366 :
367 : end subroutine cloud_particle_sedimentation_run
368 :
369 : ! Compute fluxes at interior points
370 140784 : subroutine getflx(ncol, pver, pverp, &
371 140784 : xw, phi, vel, deltat, flux, &
372 0 : errmsg, errflg)
373 : !.....xw1.......xw2.......xw3.......xw4.......xw5.......xw6
374 : !....psiw1.....psiw2.....psiw3.....psiw4.....psiw5.....psiw6
375 : !....velw1.....velw2.....velw3.....velw4.....velw5.....velw6
376 : !.........phi1......phi2.......phi3.....phi4.......phi5.......
377 :
378 : integer, intent(in) :: ncol
379 : integer, intent(in) :: pver
380 : integer, intent(in) :: pverp
381 : real(kind_phys), intent(in) :: xw(ncol, pverp)
382 : real(kind_phys), intent(in) :: phi(ncol, pver)
383 : real(kind_phys), intent(in) :: vel(ncol, pverp)
384 : real(kind_phys), intent(in) :: deltat
385 :
386 : real(kind_phys), intent(out) :: flux(ncol, pverp)
387 : character(len=512), intent(out) :: errmsg ! error message
388 : integer, intent(out) :: errflg ! error flag
389 :
390 : integer :: i, k
391 281568 : real(kind_phys) :: psi(ncol, pverp)
392 140784 : real(kind_phys) :: fdot(ncol, pverp)
393 : real(kind_phys) :: xx(ncol)
394 281568 : real(kind_phys) :: fxdot(ncol)
395 281568 : real(kind_phys) :: fxdd(ncol)
396 281568 : real(kind_phys) :: psistar(ncol)
397 281568 : real(kind_phys) :: xxk(ncol, pver)
398 :
399 2181984 : do i = 1, ncol
400 : ! integral of phi
401 2041200 : psi(i, 1) = 0._kind_phys
402 : ! fluxes at boundaries
403 2041200 : flux(i, 1) = 0._kind_phys
404 2181984 : flux(i, pverp) = 0._kind_phys
405 : end do
406 :
407 : ! integral function
408 3801168 : do k = 2, pverp
409 56872368 : do i = 1, ncol
410 56731584 : psi(i, k) = phi(i, k - 1)*(xw(i, k) - xw(i, k - 1)) + psi(i, k - 1)
411 : end do
412 : end do
413 :
414 : ! Calculate the derivatives for the interpolating polynomial
415 140784 : call cfdotmc(ncol, pver, pverp, xw, psi, fdot)
416 :
417 : ! NEW WAY
418 : ! calculate fluxes at interior pts
419 3660384 : do k = 2, pver
420 54690384 : do i = 1, ncol
421 54549600 : xxk(i, k) = xw(i, k) - vel(i, k)*deltat
422 : end do
423 : end do
424 3660384 : do k = 2, pver
425 3519600 : call cfint2(ncol, pverp, xw, psi, fdot, xxk(1, k), fxdot, fxdd, psistar, errmsg, errflg)
426 3519600 : if(errflg /= 0) then
427 0 : return
428 : endif
429 58209984 : do i = 1, ncol
430 54549600 : flux(i, k) = (psi(i, k) - psistar(i))
431 : end do
432 : end do
433 : end subroutine getflx
434 :
435 : ! Vertical level interpolation
436 3519600 : subroutine cfint2(ncol, pverp, &
437 3519600 : x, f, fdot, xin, fxdot, fxdd, psistar, &
438 0 : errmsg, errflg)
439 :
440 : integer, intent(in) :: ncol
441 : integer, intent(in) :: pverp
442 : real(kind_phys), intent(in) :: x(ncol, pverp)
443 : real(kind_phys), intent(in) :: f(ncol, pverp)
444 : real(kind_phys), intent(in) :: fdot(ncol, pverp)
445 : real(kind_phys), intent(in) :: xin(ncol)
446 :
447 : real(kind_phys), intent(out) :: fxdot(ncol)
448 : real(kind_phys), intent(out) :: fxdd(ncol)
449 : real(kind_phys), intent(out) :: psistar(ncol)
450 : character(len=512), intent(out) :: errmsg ! error message
451 : integer, intent(out) :: errflg ! error flag
452 :
453 : integer :: i, k
454 7039200 : integer :: intz(ncol)
455 : real(kind_phys) :: dx
456 : real(kind_phys) :: s
457 : real(kind_phys) :: c2
458 : real(kind_phys) :: c3
459 : real(kind_phys) :: xx
460 : real(kind_phys) :: xinf
461 : real(kind_phys) :: psi1, psi2, psi3, psim
462 : real(kind_phys) :: cfint
463 : real(kind_phys) :: cfnew
464 3519600 : real(kind_phys) :: xins(ncol)
465 : logical :: found_error
466 :
467 54549600 : do i = 1, ncol
468 51030000 : xins(i) = medan(x(i, 1), xin(i), x(i, pverp))
469 54549600 : intz(i) = 0
470 : end do
471 :
472 : ! first find the interval
473 95029200 : do k = 1, pverp - 1
474 1421809200 : do i = 1, ncol
475 1418289600 : if ((xins(i) - x(i, k))*(x(i, k + 1) - xins(i)) .ge. 0) then
476 91541169 : intz(i) = k
477 : end if
478 : end do
479 : end do
480 :
481 3519600 : found_error = .false.
482 54549600 : do i = 1, ncol
483 54549600 : if (intz(i) .eq. 0._kind_phys) found_error = .true.
484 : end do
485 3519600 : if (found_error) then
486 0 : do i = 1, ncol
487 0 : if (intz(i) .eq. 0._kind_phys) then
488 0 : write(errmsg,'(a,i0)') 'cloud_particle_sedimentation/cfint2: interval was not found for column ', i
489 0 : errflg = 1
490 : end if
491 : end do
492 : end if
493 :
494 : ! now interpolate
495 54549600 : do i = 1, ncol
496 51030000 : k = intz(i)
497 51030000 : dx = (x(i, k + 1) - x(i, k))
498 51030000 : s = (f(i, k + 1) - f(i, k))/dx
499 51030000 : c2 = (3*s - 2*fdot(i, k) - fdot(i, k + 1))/dx
500 51030000 : c3 = (fdot(i, k) + fdot(i, k + 1) - 2*s)/dx**2
501 51030000 : xx = (xins(i) - x(i, k))
502 51030000 : fxdot(i) = (3*c3*xx + 2*c2)*xx + fdot(i, k)
503 51030000 : fxdd(i) = 6*c3*xx + 2*c2
504 51030000 : cfint = ((c3*xx + c2)*xx + fdot(i, k))*xx + f(i, k)
505 :
506 : ! limit the interpolant
507 51030000 : psi1 = f(i, k) + (f(i, k + 1) - f(i, k))*xx/dx
508 51030000 : if (k .eq. 1) then
509 0 : psi2 = f(i, 1)
510 : else
511 51030000 : psi2 = f(i, k) + (f(i, k) - f(i, k - 1))*xx/(x(i, k) - x(i, k - 1))
512 : end if
513 51030000 : if (k + 1 .eq. pverp) then
514 1148592 : psi3 = f(i, pverp)
515 : else
516 49881408 : psi3 = f(i, k + 1) - (f(i, k + 2) - f(i, k + 1))*(dx - xx)/(x(i, k + 2) - x(i, k + 1))
517 : end if
518 51030000 : psim = medan(psi1, psi2, psi3)
519 51030000 : cfnew = medan(cfint, psi1, psim)
520 : !if (abs(cfnew - cfint)/(abs(cfnew) + abs(cfint) + 1.e-36_kind_phys) .gt. .03_kind_phys) then
521 : ! CHANGE THIS BACK LATER!!!
522 : ! $ .gt..1) then
523 : ! UNCOMMENT THIS LATER!!!
524 : ! write(iulog,*) ' cfint2 limiting important ', cfint, cfnew
525 : !end if
526 54549600 : psistar(i) = cfnew
527 : end do
528 3519600 : end subroutine cfint2
529 :
530 : ! Calculate the derivative for the interpolating polynomial, multi-column version
531 140784 : subroutine cfdotmc(ncol, pver, pverp, x, f, fdot)
532 : ! assumed variable distribution
533 : ! x1.......x2.......x3.......x4.......x5.......x6 1,pverp points
534 : ! f1.......f2.......f3.......f4.......f5.......f6 1,pverp points
535 : ! ...sh1.......sh2......sh3......sh4......sh5.... 1,pver points
536 : ! .........d2.......d3.......d4.......d5......... 2,pver points
537 : ! .........s2.......s3.......s4.......s5......... 2,pver points
538 : ! .............dh2......dh3......dh4............. 2,pver-1 points
539 : ! .............eh2......eh3......eh4............. 2,pver-1 points
540 : ! ..................e3.......e4.................. 3,pver-1 points
541 : ! .................ppl3......ppl4................ 3,pver-1 points
542 : ! .................ppr3......ppr4................ 3,pver-1 points
543 : ! .................t3........t4.................. 3,pver-1 points
544 : ! ................fdot3.....fdot4................ 3,pver-1 points
545 : integer, intent(in) :: ncol
546 : integer, intent(in) :: pver
547 : integer, intent(in) :: pverp
548 : real(kind_phys), intent(in) :: x(ncol, pverp)
549 : real(kind_phys), intent(in) :: f(ncol, pverp)
550 :
551 : real(kind_phys), intent(out) :: fdot(ncol, pverp) ! derivative at nodes
552 :
553 : integer :: i, k
554 : real(kind_phys) :: a, b, c ! work var
555 281568 : real(kind_phys) :: s(ncol, pverp) ! first divided differences at nodes
556 281568 : real(kind_phys) :: sh(ncol, pverp) ! first divided differences between nodes
557 281568 : real(kind_phys) :: d(ncol, pverp) ! second divided differences at nodes
558 281568 : real(kind_phys) :: dh(ncol, pverp) ! second divided differences between nodes
559 281568 : real(kind_phys) :: e(ncol, pverp) ! third divided differences at nodes
560 281568 : real(kind_phys) :: eh(ncol, pverp) ! third divided differences between nodes
561 : real(kind_phys) :: pp ! p prime
562 281568 : real(kind_phys) :: ppl(ncol, pverp) ! p prime on left
563 281568 : real(kind_phys) :: ppr(ncol, pverp) ! p prime on right
564 : real(kind_phys) :: qpl
565 : real(kind_phys) :: qpr
566 : real(kind_phys) :: ttt
567 : real(kind_phys) :: t
568 : real(kind_phys) :: tmin
569 : real(kind_phys) :: tmax
570 281568 : real(kind_phys) :: delxh(ncol, pverp)
571 :
572 3801168 : do k = 1, pver
573 : ! first divided differences between nodes
574 56731584 : do i = 1, ncol
575 53071200 : delxh(i, k) = (x(i, k + 1) - x(i, k))
576 56731584 : sh(i, k) = (f(i, k + 1) - f(i, k))/delxh(i, k)
577 : end do
578 :
579 : ! first and second divided differences at nodes
580 3801168 : if (k .ge. 2) then
581 54549600 : do i = 1, ncol
582 51030000 : d(i, k) = (sh(i, k) - sh(i, k - 1))/(x(i, k + 1) - x(i, k - 1))
583 54549600 : s(i, k) = minmod(sh(i, k), sh(i, k - 1))
584 : end do
585 : end if
586 : end do
587 :
588 : ! second and third divided diffs between nodes
589 3519600 : do k = 2, pver - 1
590 52508400 : do i = 1, ncol
591 48988800 : eh(i, k) = (d(i, k + 1) - d(i, k))/(x(i, k + 2) - x(i, k - 1))
592 52367616 : dh(i, k) = minmod(d(i, k), d(i, k + 1))
593 : end do
594 : end do
595 :
596 : ! treat the boundaries
597 2181984 : do i = 1, ncol
598 2041200 : e(i, 2) = eh(i, 2)
599 2041200 : e(i, pver) = eh(i, pver - 1)
600 : ! outside level
601 10206000 : fdot(i, 1) = sh(i, 1) - d(i, 2)*delxh(i, 1) &
602 10206000 : - eh(i, 2)*delxh(i, 1)*(x(i, 1) - x(i, 3))
603 2041200 : fdot(i, 1) = minmod(fdot(i, 1), 3*sh(i, 1))
604 16329600 : fdot(i, pverp) = sh(i, pver) + d(i, pver)*delxh(i, pver) &
605 18370800 : + eh(i, pver - 1)*delxh(i, pver)*(x(i, pverp) - x(i, pver - 1))
606 2041200 : fdot(i, pverp) = minmod(fdot(i, pverp), 3*sh(i, pver))
607 : ! one in from boundary
608 2041200 : fdot(i, 2) = sh(i, 1) + d(i, 2)*delxh(i, 1) - eh(i, 2)*delxh(i, 1)*delxh(i, 2)
609 2041200 : fdot(i, 2) = minmod(fdot(i, 2), 3*s(i, 2))
610 16329600 : fdot(i, pver) = sh(i, pver) - d(i, pver)*delxh(i, pver) &
611 18370800 : - eh(i, pver - 1)*delxh(i, pver)*delxh(i, pver - 1)
612 2181984 : fdot(i, pver) = minmod(fdot(i, pver), 3*s(i, pver))
613 : end do
614 :
615 3378816 : do k = 3, pver - 1
616 50326416 : do i = 1, ncol
617 50185632 : e(i, k) = minmod(eh(i, k), eh(i, k - 1))
618 : end do
619 : end do
620 :
621 3378816 : do k = 3, pver - 1
622 50326416 : do i = 1, ncol
623 : ! p prime at k-0.5
624 46947600 : ppl(i, k) = sh(i, k - 1) + dh(i, k - 1)*delxh(i, k - 1)
625 : ! p prime at k+0.5
626 46947600 : ppr(i, k) = sh(i, k) - dh(i, k)*delxh(i, k)
627 :
628 46947600 : t = minmod(ppl(i, k), ppr(i, k))
629 :
630 : ! derivate from parabola thru f(i,k-1), f(i,k), and f(i,k+1)
631 46947600 : pp = sh(i, k - 1) + d(i, k)*delxh(i, k - 1)
632 :
633 : ! quartic estimate of fdot
634 93895200 : fdot(i, k) = pp &
635 187790400 : - delxh(i, k - 1)*delxh(i, k) &
636 281685600 : *(eh(i, k - 1)*(x(i, k + 2) - x(i, k)) &
637 281685600 : + eh(i, k)*(x(i, k) - x(i, k - 2)) &
638 892004400 : )/(x(i, k + 2) - x(i, k - 2))
639 :
640 : ! now limit it
641 93895200 : qpl = sh(i, k - 1) &
642 469476000 : + delxh(i, k - 1)*minmod(d(i, k - 1) + e(i, k - 1)*(x(i, k) - x(i, k - 2)), &
643 610318800 : d(i, k) - e(i, k)*delxh(i, k))
644 93895200 : qpr = sh(i, k) &
645 375580800 : + delxh(i, k)*minmod(d(i, k) + e(i, k)*delxh(i, k - 1), &
646 516423600 : d(i, k + 1) + e(i, k + 1)*(x(i, k) - x(i, k + 2)))
647 :
648 46947600 : fdot(i, k) = medan(fdot(i, k), qpl, qpr)
649 :
650 46947600 : ttt = minmod(qpl, qpr)
651 46947600 : tmin = min(0._kind_phys, 3*s(i, k), 1.5_kind_phys*t, ttt)
652 46947600 : tmax = max(0._kind_phys, 3*s(i, k), 1.5_kind_phys*t, ttt)
653 :
654 50185632 : fdot(i, k) = fdot(i, k) + minmod(tmin - fdot(i, k), tmax - fdot(i, k))
655 : end do
656 : end do
657 140784 : end subroutine cfdotmc
658 :
659 589906800 : pure function minmod(a, b) result(res)
660 : real(kind_phys), intent(in) :: a, b
661 : real(kind_phys) :: res
662 589906800 : res = 0.5_kind_phys*(sign(1._kind_phys, a) + sign(1._kind_phys, b))*min(abs(a), abs(b))
663 589906800 : end function minmod
664 :
665 200037600 : pure function medan(a, b, c) result(res)
666 : real(kind_phys), intent(in) :: a, b, c
667 : real(kind_phys) :: res
668 200037600 : res = a + minmod(b - a, c - a)
669 200037600 : end function medan
670 :
671 : end module cloud_particle_sedimentation
|