[381] | 1 | MODULE dcmip2016_simple_physics_mod |
---|
| 2 | |
---|
| 3 | |
---|
| 4 | |
---|
| 5 | |
---|
| 6 | CONTAINS |
---|
| 7 | !----------------------------------------------------------------------- |
---|
| 8 | ! |
---|
| 9 | ! Date: April 26, 2016 (Version 6) |
---|
| 10 | ! |
---|
| 11 | ! Simple Physics Package |
---|
| 12 | ! |
---|
| 13 | ! SIMPLE_PHYSICS includes large-scale precipitation, surface fluxes and |
---|
| 14 | ! boundary-leyer mixing. The processes are time-split in that order. |
---|
| 15 | ! A partially implicit formulation is used to foster numerical |
---|
| 16 | ! stability. The routine assumes that the model levels are ordered |
---|
| 17 | ! in a top-down approach, e.g. level 1 denotes the uppermost full model |
---|
| 18 | ! level. |
---|
| 19 | ! |
---|
| 20 | ! This routine is based on an implementation which was developed for |
---|
| 21 | ! the NCAR Community Atmosphere Model (CAM). Adjustments for other |
---|
| 22 | ! models may be necessary. |
---|
| 23 | ! |
---|
| 24 | ! The routine provides both updates of the state variables u, v, T, q |
---|
| 25 | ! (these are local copies of u,v,T,q within this physics routine) and |
---|
| 26 | ! also collects their time tendencies. The latter might be used to |
---|
| 27 | ! couple the physics and dynamics in a process-split way. For a |
---|
| 28 | ! time-split coupling, the final state should be given to the |
---|
| 29 | ! dynamical core for the next time step. |
---|
| 30 | ! |
---|
| 31 | ! Test: 0 = Reed and Jablonowski (2011) tropical cyclone test |
---|
| 32 | ! 1 = Moist baroclinic instability test |
---|
| 33 | ! RJ2012_precip: |
---|
| 34 | ! true = Turn on Reed and Jablonowski (2012) precip scheme |
---|
| 35 | ! false = Turn off Reed and Jablonowski (2012) precip scheme |
---|
| 36 | ! TC_PBL_mod: |
---|
| 37 | ! true = Turn on George Bryan PBL mod for tropical cyclone test |
---|
| 38 | ! false = Turn off George Bryan PBL mod (i.e., run as in Reed and Jablonowski (2012)) |
---|
| 39 | ! |
---|
| 40 | ! SUBROUTINE SIMPLE_PHYSICS(pcols, pver, dtime, lat, t, q, u, v, pmid, pint, pdel, rpdel, ps, precl, test, RJ2012_precip, TC_PBL_mod) |
---|
| 41 | ! |
---|
| 42 | ! Input variables: |
---|
| 43 | ! pcols - number of atmospheric columns (#) |
---|
| 44 | ! pver - number of model levels (#) |
---|
| 45 | ! dtime - time step (s) |
---|
| 46 | ! lat - latitude (radians) |
---|
| 47 | ! t - temperature at model levels (K) |
---|
| 48 | ! q - specific humidity at model levels (gm/gm) |
---|
| 49 | ! u - zonal wind at model levels (m/s) |
---|
| 50 | ! v - meridional wind at model levels (m/s) |
---|
| 51 | ! pmid - pressure at model levels (Pa) |
---|
| 52 | ! pint - pressure at interfaces (Pa) |
---|
| 53 | ! pdel - layer thickness (Pa) |
---|
| 54 | ! rpdel - reciprocal of layer thickness (1/Pa) |
---|
| 55 | ! ps - surface pressure (Pa) |
---|
| 56 | ! test - test case to use for sea-surface temperatures |
---|
| 57 | ! RJ2012_precip - RJ2012 precip flag |
---|
| 58 | ! TC_PBL_mod - PCL modification for TC test |
---|
| 59 | ! |
---|
| 60 | ! Output variables: |
---|
| 61 | ! Increments are added into t, q, u, v, pmid, pint, pdel, rpdel and ps |
---|
| 62 | ! which are returned to the routine from which SIMPLE_PHYSICS was |
---|
| 63 | ! called. Precpitation is returned via precl. |
---|
| 64 | ! |
---|
| 65 | ! Change log: |
---|
| 66 | ! v2: removal of some NCAR CAM-specific 'use' associations |
---|
| 67 | ! v3: corrected precl(i) computation, the precipitation rate is now |
---|
| 68 | ! computed via a vertical integral, the previous single-level |
---|
| 69 | ! computation in v2 was a bug |
---|
| 70 | ! v3: corrected dtdt(i,1) computation, the term '-(i,1)' was missing |
---|
| 71 | ! the temperature variable: '-t(i,1)' |
---|
| 72 | ! v4: modified and enhanced parameter list to make the routine truly |
---|
| 73 | ! standalone, the number of columns and vertical levels have been |
---|
| 74 | ! added: pcols, pver |
---|
| 75 | ! v4: 'ncol' has been removed, 'pcols' is used instead |
---|
| 76 | ! v5: the sea surface temperature (SST) field Tsurf is now an array, |
---|
| 77 | ! the SST now depends on the latitude |
---|
| 78 | ! v5: addition of the latitude array 'lat' and the flag 'test' in the |
---|
| 79 | ! parameter list |
---|
| 80 | ! if test = 0: constant SST is used, correct setting for the |
---|
| 81 | ! tropical cyclone test |
---|
| 82 | ! if test = 1: newly added latitude-dependent SST is used, |
---|
| 83 | ! correct setting for the moist baroclinic wave test |
---|
| 84 | ! with simple-physics |
---|
| 85 | ! v6: addition of flags for a modified PBL for the TC test and |
---|
| 86 | ! to turn off large-scale condensation scheme when using Kessler physics |
---|
| 87 | ! Included virtual temperature in density calculation in PBL scheme |
---|
| 88 | ! Also, included the virtual temperature, instead of temperature, for |
---|
| 89 | ! the calculation of rho in the PBL scheme |
---|
| 90 | ! (v6_1) Minor specification and generalization fixes. |
---|
| 91 | ! |
---|
| 92 | ! Reference: Reed, K. A. and C. Jablonowski (2012), Idealized tropical cyclone |
---|
| 93 | ! simulations of intermediate complexity: A test case for AGCMs, |
---|
| 94 | ! J. Adv. Model. Earth Syst., Vol. 4, M04001, doi:10.1029/2011MS000099 |
---|
| 95 | !----------------------------------------------------------------------- |
---|
| 96 | |
---|
| 97 | SUBROUTINE SIMPLE_PHYSICS(pcols, pver, dtime, lat, t, q, u, v, pmid, pint, pdel, rpdel, ps, precl, test, RJ2012_precip, TC_PBL_mod) |
---|
| 98 | |
---|
| 99 | ! use physics_types , only: physics_dme_adjust ! This is for CESM/CAM |
---|
| 100 | ! use cam_diagnostics, only: diag_phys_writeout ! This is for CESM/CAM |
---|
| 101 | |
---|
| 102 | implicit none |
---|
| 103 | |
---|
| 104 | integer, parameter :: r8 = selected_real_kind(12) |
---|
| 105 | |
---|
| 106 | ! |
---|
| 107 | ! Input arguments - MODEL DEPENDENT |
---|
| 108 | ! |
---|
| 109 | integer, intent(in) :: pcols ! Set number of atmospheric columns |
---|
| 110 | integer, intent(in) :: pver ! Set number of model levels |
---|
| 111 | real(r8), intent(in) :: dtime ! Set model physics timestep |
---|
| 112 | real(r8), intent(in) :: lat(pcols) ! Latitude |
---|
| 113 | integer, intent(in) :: test ! Test number |
---|
| 114 | logical, intent(in) :: RJ2012_precip |
---|
| 115 | logical, intent(in) :: TC_PBL_mod |
---|
| 116 | |
---|
| 117 | ! |
---|
| 118 | ! Input/Output arguments |
---|
| 119 | ! |
---|
| 120 | ! pcols is the maximum number of vertical columns per 'chunk' of atmosphere |
---|
| 121 | ! |
---|
| 122 | real(r8), intent(inout) :: t(pcols,pver) ! Temperature at full-model level (K) |
---|
| 123 | real(r8), intent(inout) :: q(pcols,pver) ! Specific Humidity at full-model level (kg/kg) |
---|
| 124 | real(r8), intent(inout) :: u(pcols,pver) ! Zonal wind at full-model level (m/s) |
---|
| 125 | real(r8), intent(inout) :: v(pcols,pver) ! Meridional wind at full-model level (m/s) |
---|
| 126 | real(r8), intent(inout) :: pmid(pcols,pver) ! Pressure is full-model level (Pa) |
---|
| 127 | real(r8), intent(inout) :: pint(pcols,pver+1) ! Pressure at model interfaces (Pa) |
---|
| 128 | real(r8), intent(inout) :: pdel(pcols,pver) ! Layer thickness (Pa) |
---|
| 129 | real(r8), intent(in ) :: rpdel(pcols,pver) ! Reciprocal of layer thickness (1/Pa) |
---|
| 130 | real(r8), intent(inout) :: ps(pcols) ! Surface Pressue (Pa) |
---|
| 131 | |
---|
| 132 | ! |
---|
| 133 | ! Output arguments |
---|
| 134 | ! |
---|
| 135 | real(r8), intent(out) :: precl(pcols) ! Precipitation rate (m_water / s) |
---|
| 136 | |
---|
| 137 | ! |
---|
| 138 | !---------------------------Local workspace----------------------------- |
---|
| 139 | ! |
---|
| 140 | |
---|
| 141 | ! Integers for loops |
---|
| 142 | |
---|
| 143 | integer i,k ! Longitude, level indices |
---|
| 144 | |
---|
| 145 | ! Physical Constants - Many of these may be model dependent |
---|
| 146 | |
---|
| 147 | real(r8) gravit ! Gravity |
---|
| 148 | real(r8) rair ! Gas constant for dry air |
---|
| 149 | real(r8) cpair ! Specific heat of dry air |
---|
| 150 | real(r8) latvap ! Latent heat of vaporization |
---|
| 151 | real(r8) rh2o ! Gas constant for water vapor |
---|
| 152 | real(r8) epsilo ! Ratio of gas constant for dry air to that for vapor |
---|
| 153 | real(r8) zvir ! Constant for virtual temp. calc. =(rh2o/rair) - 1 |
---|
| 154 | real(r8) a ! Reference Earth's Radius (m) |
---|
| 155 | real(r8) omega ! Reference rotation rate of the Earth (s^-1) |
---|
| 156 | real(r8) pi ! pi |
---|
| 157 | |
---|
| 158 | ! Simple Physics Specific Constants |
---|
| 159 | |
---|
| 160 | !++++++++ |
---|
| 161 | real(r8) Tsurf(pcols) ! Sea Surface Temperature (constant for tropical cyclone) |
---|
| 162 | !++++++++ Tsurf needs to be dependent on latitude for the |
---|
| 163 | ! moist baroclinic wave test, adjust |
---|
| 164 | |
---|
| 165 | real(r8) SST_TC ! Sea Surface Temperature for tropical cyclone test |
---|
| 166 | real(r8) T0 ! Control temp for calculation of qsat |
---|
| 167 | real(r8) e0 ! Saturation vapor pressure at T0 for calculation of qsat |
---|
| 168 | real(r8) rhow ! Density of Liquid Water |
---|
| 169 | real(r8) p0 ! Constant for calculation of potential temperature |
---|
| 170 | real(r8) Cd0 ! Constant for calculating Cd from Smith and Vogl 2008 |
---|
| 171 | real(r8) Cd1 ! Constant for calculating Cd from Smith and Vogl 2008 |
---|
| 172 | real(r8) Cm ! Constant for calculating Cd from Smith and Vogl 2008 |
---|
| 173 | real(r8) v20 ! Threshold wind speed for calculating Cd from Smith and Vogl 2008 |
---|
| 174 | real(r8) C ! Drag coefficient for sensible heat and evaporation |
---|
| 175 | real(r8) T00 ! Horizontal mean T at surface for moist baro test |
---|
| 176 | real(r8) u0 ! Zonal wind constant for moist baro test |
---|
| 177 | real(r8) latw ! halfwidth for for baro test |
---|
| 178 | real(r8) eta0 ! Center of jets (hybrid) for baro test |
---|
| 179 | real(r8) etav ! Auxiliary variable for baro test |
---|
| 180 | real(r8) q0 ! Maximum specific humidity for baro test |
---|
| 181 | real(r8) kappa ! von Karman constant |
---|
| 182 | |
---|
| 183 | ! Physics Tendency Arrays |
---|
| 184 | real(r8) dtdt(pcols,pver) ! Temperature tendency |
---|
| 185 | real(r8) dqdt(pcols,pver) ! Specific humidity tendency |
---|
| 186 | real(r8) dudt(pcols,pver) ! Zonal wind tendency |
---|
| 187 | real(r8) dvdt(pcols,pver) ! Meridional wind tendency |
---|
| 188 | |
---|
| 189 | ! Temporary variables for tendency calculations |
---|
| 190 | |
---|
| 191 | real(r8) tmp ! Temporary |
---|
| 192 | real(r8) qsat ! Saturation vapor pressure |
---|
| 193 | real(r8) qsats ! Saturation vapor pressure of SST |
---|
| 194 | |
---|
| 195 | ! Variables for Boundary Layer Calculation |
---|
| 196 | |
---|
| 197 | real(r8) wind(pcols) ! Magnitude of Wind |
---|
| 198 | real(r8) Cd(pcols) ! Drag coefficient for momentum |
---|
| 199 | real(r8) Km(pcols,pver+1) ! Eddy diffusivity for boundary layer calculations |
---|
| 200 | real(r8) Ke(pcols,pver+1) ! Eddy diffusivity for boundary layer calculations |
---|
| 201 | real(r8) rho ! Density at lower/upper interface |
---|
| 202 | real(r8) za(pcols) ! Heights at midpoints of first model level |
---|
| 203 | real(r8) zi(pcols,pver+1) ! Heights at model interfaces |
---|
| 204 | real(r8) dlnpint ! Used for calculation of heights |
---|
| 205 | real(r8) pbltop ! Top of boundary layer |
---|
| 206 | real(r8) zpbltop ! Top of boundary layer for George Bryan Modifcation |
---|
| 207 | real(r8) pblconst ! Constant for the calculation of the decay of diffusivity |
---|
| 208 | real(r8) CA(pcols,pver) ! Matrix Coefficents for PBL Scheme |
---|
| 209 | real(r8) CC(pcols,pver) ! Matrix Coefficents for PBL Scheme |
---|
| 210 | real(r8) CE(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 211 | real(r8) CAm(pcols,pver) ! Matrix Coefficents for PBL Scheme |
---|
| 212 | real(r8) CCm(pcols,pver) ! Matrix Coefficents for PBL Scheme |
---|
| 213 | real(r8) CEm(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 214 | real(r8) CFu(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 215 | real(r8) CFv(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 216 | real(r8) CFt(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 217 | real(r8) CFq(pcols,pver+1) ! Matrix Coefficents for PBL Scheme |
---|
| 218 | |
---|
| 219 | |
---|
| 220 | ! Variable for Dry Mass Adjustment, this dry air adjustment is necessary to |
---|
| 221 | ! conserve the mass of the dry air |
---|
| 222 | |
---|
| 223 | real(r8) qini(pcols,pver) ! Initial specific humidity |
---|
| 224 | |
---|
| 225 | !=============================================================================== |
---|
| 226 | ! |
---|
| 227 | ! Physical Constants - MAY BE MODEL DEPENDENT |
---|
| 228 | ! |
---|
| 229 | !=============================================================================== |
---|
| 230 | gravit = 9.80616_r8 ! Gravity (9.80616 m/s^2) |
---|
| 231 | rair = 287.0_r8 ! Gas constant for dry air: 287 J/(kg K) |
---|
| 232 | cpair = 1.0045e3_r8 ! Specific heat of dry air: here we use 1004.5 J/(kg K) |
---|
| 233 | latvap = 2.5e6_r8 ! Latent heat of vaporization (J/kg) |
---|
| 234 | rh2o = 461.5_r8 ! Gas constant for water vapor: 461.5 J/(kg K) |
---|
| 235 | epsilo = rair/rh2o ! Ratio of gas constant for dry air to that for vapor |
---|
| 236 | zvir = (rh2o/rair) - 1._r8 ! Constant for virtual temp. calc. =(rh2o/rair) - 1 is approx. 0.608 |
---|
| 237 | a = 6371220.0_r8 ! Reference Earth's Radius (m) |
---|
| 238 | omega = 7.29212d-5 ! Reference rotation rate of the Earth (s^-1) |
---|
| 239 | pi = 4._r8*atan(1._r8) ! pi |
---|
| 240 | |
---|
| 241 | !=============================================================================== |
---|
| 242 | ! |
---|
| 243 | ! Local Constants for Simple Physics |
---|
| 244 | ! |
---|
| 245 | !=============================================================================== |
---|
| 246 | C = 0.0011_r8 ! From Simth and Vogl 2008 |
---|
| 247 | SST_TC = 302.15_r8 ! Constant Value for SST |
---|
| 248 | T0 = 273.16_r8 ! control temp for calculation of qsat |
---|
| 249 | e0 = 610.78_r8 ! saturation vapor pressure at T0 for calculation of qsat |
---|
| 250 | rhow = 1000.0_r8 ! Density of Liquid Water |
---|
| 251 | Cd0 = 0.0007_r8 ! Constant for Cd calc. Simth and Vogl 2008 |
---|
| 252 | Cd1 = 0.000065_r8 ! Constant for Cd calc. Simth and Vogl 2008 |
---|
| 253 | Cm = 0.002_r8 ! Constant for Cd calc. Simth and Vogl 2008 |
---|
| 254 | v20 = 20.0_r8 ! Threshold wind speed for calculating Cd from Smith and Vogl 2008 |
---|
| 255 | p0 = 100000.0_r8 ! Constant for potential temp calculation |
---|
| 256 | pbltop = 85000._r8 ! Top of boundary layer in p |
---|
| 257 | zpbltop = 1000._r8 ! Top of boundary layer in z |
---|
| 258 | pblconst = 10000._r8 ! Constant for the calculation of the decay of diffusivity |
---|
| 259 | T00 = 288.0_r8 ! Horizontal mean T at surface for moist baro test |
---|
| 260 | u0 = 35.0_r8 ! Zonal wind constant for moist baro test |
---|
| 261 | latw = 2.0_r8*pi/9.0_r8 ! Halfwidth for for baro test |
---|
| 262 | eta0 = 0.252_r8 ! Center of jets (hybrid) for baro test |
---|
| 263 | etav = (1._r8-eta0)*0.5_r8*pi ! Auxiliary variable for baro test |
---|
| 264 | q0 = 0.021_r8 ! Maximum specific humidity for baro test |
---|
| 265 | kappa = 0.4_r8 ! von Karman constant |
---|
| 266 | |
---|
| 267 | !=============================================================================== |
---|
| 268 | ! |
---|
| 269 | ! Definition of local arrays |
---|
| 270 | ! |
---|
| 271 | !=============================================================================== |
---|
| 272 | ! |
---|
| 273 | ! Calculate hydrostatic height |
---|
| 274 | ! |
---|
| 275 | do i=1,pcols |
---|
| 276 | dlnpint = log(ps(i)) - log(pint(i,pver)) ! ps(i) is identical to pint(i,pver+1), note: this is the correct sign (corrects typo in JAMES paper) |
---|
| 277 | za(i) = rair/gravit*t(i,pver)*(1._r8+zvir*q(i,pver))*0.5_r8*dlnpint |
---|
| 278 | zi(i,pver+1) = 0.0_r8 |
---|
| 279 | end do |
---|
| 280 | ! |
---|
| 281 | ! Set Initial Specific Humidity |
---|
| 282 | ! |
---|
| 283 | qini(:pcols,:pver) = q(:pcols,:pver) |
---|
| 284 | ! |
---|
| 285 | ! Set Sea Surface Temperature (constant for tropical cyclone) |
---|
| 286 | ! Tsurf needs to be dependent on latitude for moist baroclinic wave test |
---|
| 287 | ! Tsurf needs to be constant for tropical cyclone test |
---|
| 288 | ! |
---|
| 289 | if (test .eq. 1) then ! Moist Baroclinic Wave Test |
---|
| 290 | do i=1,pcols |
---|
| 291 | Tsurf(i) = (T00 + pi*u0/rair * 1.5_r8 * sin(etav) * (cos(etav))**0.5_r8 * & |
---|
| 292 | ((-2._r8*(sin(lat(i)))**6 * ((cos(lat(i)))**2 + 1._r8/3._r8) + 10._r8/63._r8)* & |
---|
| 293 | u0 * (cos(etav))**1.5_r8 + & |
---|
| 294 | (8._r8/5._r8*(cos(lat(i)))**3 * ((sin(lat(i)))**2 + 2._r8/3._r8) - pi/4._r8)*a*omega*0.5_r8 ))/ & |
---|
| 295 | (1._r8+zvir*q0*exp(-(lat(i)/latw)**4)) |
---|
| 296 | |
---|
| 297 | end do |
---|
| 298 | else if (test .eq. 0) then ! Tropical Cyclone Test |
---|
| 299 | do i=1,pcols |
---|
| 300 | Tsurf(i) = SST_TC |
---|
| 301 | end do |
---|
| 302 | end if |
---|
| 303 | |
---|
| 304 | !=============================================================================== |
---|
| 305 | ! |
---|
| 306 | ! Set initial physics time tendencies and precipitation field to zero |
---|
| 307 | ! |
---|
| 308 | !=============================================================================== |
---|
| 309 | dtdt(:pcols,:pver) = 0._r8 ! initialize temperature tendency with zero |
---|
| 310 | dqdt(:pcols,:pver) = 0._r8 ! initialize specific humidity tendency with zero |
---|
| 311 | dudt(:pcols,:pver) = 0._r8 ! initialize zonal wind tendency with zero |
---|
| 312 | dvdt(:pcols,:pver) = 0._r8 ! initialize meridional wind tendency with zero |
---|
| 313 | precl(:pcols) = 0._r8 ! initialize precipitation rate with zero |
---|
| 314 | |
---|
| 315 | !=============================================================================== |
---|
| 316 | ! |
---|
| 317 | ! Large-Scale Condensation and Precipitation |
---|
| 318 | ! |
---|
| 319 | !=============================================================================== |
---|
| 320 | |
---|
| 321 | if (RJ2012_precip) then |
---|
| 322 | ! |
---|
| 323 | ! Calculate Tendencies |
---|
| 324 | ! |
---|
| 325 | do k=1,pver |
---|
| 326 | do i=1,pcols |
---|
| 327 | qsat = epsilo*e0/pmid(i,k)*exp(-latvap/rh2o*((1._r8/t(i,k))-1._r8/T0)) |
---|
| 328 | if (q(i,k) > qsat) then |
---|
| 329 | tmp = 1._r8/dtime*(q(i,k)-qsat)/(1._r8+(latvap/cpair)*(epsilo*latvap*qsat/(rair*t(i,k)**2))) |
---|
| 330 | dtdt(i,k) = dtdt(i,k)+latvap/cpair*tmp |
---|
| 331 | dqdt(i,k) = dqdt(i,k)-tmp |
---|
| 332 | precl(i) = precl(i)+tmp*pdel(i,k)/(gravit*rhow) |
---|
| 333 | end if |
---|
| 334 | end do |
---|
| 335 | end do |
---|
| 336 | ! |
---|
| 337 | ! Update moisture and temperature fields from Larger-Scale Precipitation Scheme |
---|
| 338 | ! |
---|
| 339 | do k=1,pver |
---|
| 340 | do i=1,pcols |
---|
| 341 | t(i,k) = t(i,k) + dtdt(i,k)*dtime |
---|
| 342 | q(i,k) = q(i,k) + dqdt(i,k)*dtime |
---|
| 343 | end do |
---|
| 344 | end do |
---|
| 345 | |
---|
| 346 | !=============================================================================== |
---|
| 347 | ! Send variables to history file - THIS PROCESS WILL BE MODEL SPECIFIC |
---|
| 348 | ! |
---|
| 349 | ! note: The variables, as done in many GCMs, are written to the history file |
---|
| 350 | ! after the moist physics process. This ensures that the moisture fields |
---|
| 351 | ! are somewhat in equilibrium. |
---|
| 352 | !=============================================================================== |
---|
| 353 | ! call diag_phys_writeout(state) ! This is for CESM/CAM |
---|
| 354 | |
---|
| 355 | end if |
---|
| 356 | |
---|
| 357 | !=============================================================================== |
---|
| 358 | ! |
---|
| 359 | ! Turbulent mixing coefficients for the PBL mixing of horizontal momentum, |
---|
| 360 | ! sensible heat and latent heat |
---|
| 361 | ! |
---|
| 362 | ! We are using Simplified Ekman theory to compute the diffusion coefficients |
---|
| 363 | ! Kx for the boundary-layer mixing. The Kx values are calculated at each time step |
---|
| 364 | ! and in each column. |
---|
| 365 | ! |
---|
| 366 | !===============================================================================! |
---|
| 367 | ! Compute magnitude of the wind and drag coeffcients for turbulence scheme |
---|
| 368 | ! |
---|
| 369 | do i=1,pcols |
---|
| 370 | wind(i) = sqrt(u(i,pver)**2+v(i,pver)**2) |
---|
| 371 | end do |
---|
| 372 | do i=1,pcols |
---|
| 373 | if( wind(i) .lt. v20) then |
---|
| 374 | Cd(i) = Cd0+Cd1*wind(i) |
---|
| 375 | else |
---|
| 376 | Cd(i) = Cm |
---|
| 377 | endif |
---|
| 378 | end do |
---|
| 379 | |
---|
| 380 | if (TC_PBL_mod) then !Bryan TC PBL Modification |
---|
| 381 | do k=pver,1,-1 |
---|
| 382 | do i=1,pcols |
---|
| 383 | dlnpint = log(pint(i,k+1)) - log(pint(i,k)) |
---|
| 384 | zi(i,k) = zi(i,k+1)+rair/gravit*t(i,k)*(1._r8+zvir*q(i,k))*dlnpint |
---|
| 385 | if( zi(i,k) .le. zpbltop) then |
---|
| 386 | Km(i,k) = kappa*sqrt(Cd(i))*wind(i)*zi(i,k)*(1._r8-zi(i,k)/zpbltop)*(1._r8-zi(i,k)/zpbltop) |
---|
| 387 | Ke(i,k) = kappa*sqrt(C)*wind(i)*zi(i,k)*(1._r8-zi(i,k)/zpbltop)*(1._r8-zi(i,k)/zpbltop) |
---|
| 388 | else |
---|
| 389 | Km(i,k) = 0.0_r8 |
---|
| 390 | Ke(i,k) = 0.0_r8 |
---|
| 391 | end if |
---|
| 392 | end do |
---|
| 393 | end do |
---|
| 394 | else ! Reed and Jablonowski (2012) Configuration |
---|
| 395 | do k=1,pver |
---|
| 396 | do i=1,pcols |
---|
| 397 | if( pint(i,k) .ge. pbltop) then |
---|
| 398 | Km(i,k) = Cd(i)*wind(i)*za(i) |
---|
| 399 | Ke(i,k) = C*wind(i)*za(i) |
---|
| 400 | else |
---|
| 401 | Km(i,k) = Cd(i)*wind(i)*za(i)*exp(-(pbltop-pint(i,k))**2/(pblconst)**2) |
---|
| 402 | Ke(i,k) = C*wind(i)*za(i)*exp(-(pbltop-pint(i,k))**2/(pblconst)**2) |
---|
| 403 | end if |
---|
| 404 | end do |
---|
| 405 | end do |
---|
| 406 | end if |
---|
| 407 | !=============================================================================== |
---|
| 408 | ! Update the state variables u, v, t, q with the surface fluxes at the |
---|
| 409 | ! lowest model level, this is done with an implicit approach |
---|
| 410 | ! see Reed and Jablonowski (JAMES, 2012) |
---|
| 411 | ! |
---|
| 412 | ! Sea Surface Temperature Tsurf is constant for tropical cyclone test 5-1 |
---|
| 413 | ! Tsurf needs to be dependent on latitude for the |
---|
| 414 | ! moist baroclinic wave test |
---|
| 415 | !=============================================================================== |
---|
| 416 | do i=1,pcols |
---|
| 417 | qsats = epsilo*e0/ps(i)*exp(-latvap/rh2o*((1._r8/Tsurf(i))-1._r8/T0)) |
---|
| 418 | dudt(i,pver) = dudt(i,pver) + (u(i,pver)/(1._r8+Cd(i)*wind(i)*dtime/za(i))-u(i,pver))/dtime |
---|
| 419 | dvdt(i,pver) = dvdt(i,pver) + (v(i,pver)/(1._r8+Cd(i)*wind(i)*dtime/za(i))-v(i,pver))/dtime |
---|
| 420 | u(i,pver) = u(i,pver)/(1._r8+Cd(i)*wind(i)*dtime/za(i)) |
---|
| 421 | v(i,pver) = v(i,pver)/(1._r8+Cd(i)*wind(i)*dtime/za(i)) |
---|
| 422 | dtdt(i,pver) = dtdt(i,pver) +((t(i,pver)+C*wind(i)*Tsurf(i)*dtime/za(i))/(1._r8+C*wind(i)*dtime/za(i))-t(i,pver))/dtime |
---|
| 423 | t(i,pver) = (t(i,pver)+C*wind(i)*Tsurf(i)*dtime/za(i))/(1._r8+C*wind(i)*dtime/za(i)) |
---|
| 424 | dqdt(i,pver) = dqdt(i,pver) +((q(i,pver)+C*wind(i)*qsats*dtime/za(i))/(1._r8+C*wind(i)*dtime/za(i))-q(i,pver))/dtime |
---|
| 425 | q(i,pver) = (q(i,pver)+C*wind(i)*qsats*dtime/za(i))/(1._r8+C*wind(i)*dtime/za(i)) |
---|
| 426 | end do |
---|
| 427 | |
---|
| 428 | !=============================================================================== |
---|
| 429 | ! Boundary layer mixing, see Reed and Jablonowski (JAMES, 2012) |
---|
| 430 | !=============================================================================== |
---|
| 431 | ! Calculate Diagonal Variables for Implicit PBL Scheme |
---|
| 432 | ! |
---|
| 433 | do k=1,pver-1 |
---|
| 434 | do i=1,pcols |
---|
| 435 | rho = (pint(i,k+1)/(rair*(t(i,k+1)*(1._r8+zvir*q(i,k+1))+t(i,k)*(1._r8+zvir*q(i,k)))/2.0_r8)) |
---|
| 436 | CAm(i,k) = rpdel(i,k)*dtime*gravit*gravit*Km(i,k+1)*rho*rho/(pmid(i,k+1)-pmid(i,k)) |
---|
| 437 | CCm(i,k+1) = rpdel(i,k+1)*dtime*gravit*gravit*Km(i,k+1)*rho*rho/(pmid(i,k+1)-pmid(i,k)) |
---|
| 438 | CA(i,k) = rpdel(i,k)*dtime*gravit*gravit*Ke(i,k+1)*rho*rho/(pmid(i,k+1)-pmid(i,k)) |
---|
| 439 | CC(i,k+1) = rpdel(i,k+1)*dtime*gravit*gravit*Ke(i,k+1)*rho*rho/(pmid(i,k+1)-pmid(i,k)) |
---|
| 440 | end do |
---|
| 441 | end do |
---|
| 442 | do i=1,pcols |
---|
| 443 | CAm(i,pver) = 0._r8 |
---|
| 444 | CCm(i,1) = 0._r8 |
---|
| 445 | CEm(i,pver+1) = 0._r8 |
---|
| 446 | CA(i,pver) = 0._r8 |
---|
| 447 | CC(i,1) = 0._r8 |
---|
| 448 | CE(i,pver+1) = 0._r8 |
---|
| 449 | CFu(i,pver+1) = 0._r8 |
---|
| 450 | CFv(i,pver+1) = 0._r8 |
---|
| 451 | CFt(i,pver+1) = 0._r8 |
---|
| 452 | CFq(i,pver+1) = 0._r8 |
---|
| 453 | end do |
---|
| 454 | do i=1,pcols |
---|
| 455 | do k=pver,1,-1 |
---|
| 456 | CE(i,k) = CC(i,k)/(1._r8+CA(i,k)+CC(i,k)-CA(i,k)*CE(i,k+1)) |
---|
| 457 | CEm(i,k) = CCm(i,k)/(1._r8+CAm(i,k)+CCm(i,k)-CAm(i,k)*CEm(i,k+1)) |
---|
| 458 | CFu(i,k) = (u(i,k)+CAm(i,k)*CFu(i,k+1))/(1._r8+CAm(i,k)+CCm(i,k)-CAm(i,k)*CEm(i,k+1)) |
---|
| 459 | CFv(i,k) = (v(i,k)+CAm(i,k)*CFv(i,k+1))/(1._r8+CAm(i,k)+CCm(i,k)-CAm(i,k)*CEm(i,k+1)) |
---|
| 460 | CFt(i,k) = ((p0/pmid(i,k))**(rair/cpair)*t(i,k)+CA(i,k)*CFt(i,k+1))/(1._r8+CA(i,k)+CC(i,k)-CA(i,k)*CE(i,k+1)) |
---|
| 461 | CFq(i,k) = (q(i,k)+CA(i,k)*CFq(i,k+1))/(1._r8+CA(i,k)+CC(i,k)-CA(i,k)*CE(i,k+1)) |
---|
| 462 | end do |
---|
| 463 | end do |
---|
| 464 | ! |
---|
| 465 | ! Calculate the updated temperaure and specific humidity and wind tendencies |
---|
| 466 | ! |
---|
| 467 | ! First we need to calculate the tendencies at the top model level |
---|
| 468 | ! |
---|
| 469 | do i=1,pcols |
---|
| 470 | dudt(i,1) = dudt(i,1)+(CFu(i,1)-u(i,1))/dtime |
---|
| 471 | dvdt(i,1) = dvdt(i,1)+(CFv(i,1)-v(i,1))/dtime |
---|
| 472 | u(i,1) = CFu(i,1) |
---|
| 473 | v(i,1) = CFv(i,1) |
---|
| 474 | dtdt(i,1) = dtdt(i,1)+(CFt(i,1)*(pmid(i,1)/p0)**(rair/cpair)-t(i,1))/dtime |
---|
| 475 | t(i,1) = CFt(i,1)*(pmid(i,1)/p0)**(rair/cpair) |
---|
| 476 | dqdt(i,1) = dqdt(i,1)+(CFq(i,1)-q(i,1))/dtime |
---|
| 477 | q(i,1) = CFq(i,1) |
---|
| 478 | end do |
---|
| 479 | |
---|
| 480 | do i=1,pcols |
---|
| 481 | do k=2,pver |
---|
| 482 | dudt(i,k) = dudt(i,k)+(CEm(i,k)*u(i,k-1)+CFu(i,k)-u(i,k))/dtime |
---|
| 483 | dvdt(i,k) = dvdt(i,k)+(CEm(i,k)*v(i,k-1)+CFv(i,k)-v(i,k))/dtime |
---|
| 484 | u(i,k) = CEm(i,k)*u(i,k-1)+CFu(i,k) |
---|
| 485 | v(i,k) = CEm(i,k)*v(i,k-1)+CFv(i,k) |
---|
[402] | 486 | dtdt(i,k) = dtdt(i,k)+ & |
---|
| 487 | ((CE(i,k)*t(i,k-1)*(p0/pmid(i,k-1))**(rair/cpair)+CFt(i,k))*(pmid(i,k)/p0)**(rair/cpair)-t(i,k))/dtime |
---|
[381] | 488 | t(i,k) = (CE(i,k)*t(i,k-1)*(p0/pmid(i,k-1))**(rair/cpair)+CFt(i,k))*(pmid(i,k)/p0)**(rair/cpair) |
---|
| 489 | dqdt(i,k) = dqdt(i,k)+(CE(i,k)*q(i,k-1)+CFq(i,k)-q(i,k))/dtime |
---|
| 490 | q(i,k) = CE(i,k)*q(i,k-1)+CFq(i,k) |
---|
| 491 | end do |
---|
| 492 | end do |
---|
| 493 | |
---|
| 494 | !=============================================================================== |
---|
| 495 | ! |
---|
| 496 | ! Dry Mass Adjustment - THIS PROCESS WILL BE MODEL SPECIFIC |
---|
| 497 | ! |
---|
| 498 | ! note: Care needs to be taken to ensure that the model conserves the total |
---|
| 499 | ! dry air mass. Add your own routine here. |
---|
| 500 | !=============================================================================== |
---|
| 501 | ! call physics_dme_adjust(state, tend, qini, dtime) ! This is for CESM/CAM |
---|
| 502 | |
---|
| 503 | return |
---|
| 504 | end subroutine SIMPLE_PHYSICS |
---|
| 505 | |
---|
| 506 | |
---|
| 507 | END MODULE dcmip2016_simple_physics_mod |
---|