[1033] | 1 | !======================================================================= |
---|
[1032] | 2 | ! CALLTHERM_INTERFACE |
---|
[1033] | 3 | !======================================================================= |
---|
[1032] | 4 | ! Main interface to the Martian thermal plume model |
---|
| 5 | ! This interface handles sub-timesteps for this model |
---|
| 6 | ! A call to this interface must be inserted in the main 'physics' routine |
---|
| 7 | ! NB: for information: |
---|
| 8 | ! In the Mars LMD-GCM, the thermal plume model is called after |
---|
| 9 | ! the vertical turbulent mixing scheme (Mellor and Yamada) |
---|
| 10 | ! and the surface layer scheme (Richardson-based surface layer + subgrid gustiness) |
---|
| 11 | ! Other routines called before the thermals model are |
---|
| 12 | ! radiative transfer and (orographic) gravity wave drag. |
---|
| 13 | ! ----------------------------------------------------------------------- |
---|
[1033] | 14 | ! Author : A. Colaitis 2011-01-05 (with updates 2011-2013) |
---|
| 15 | ! after C. Rio and F. Hourdin |
---|
| 16 | ! Institution : Laboratoire de Meteorologie Dynamique (LMD) Paris, France |
---|
| 17 | ! ----------------------------------------------------------------------- |
---|
| 18 | ! Corresponding author : A. Spiga aymeric.spiga_AT_upmc.fr |
---|
| 19 | ! ----------------------------------------------------------------------- |
---|
[1032] | 20 | ! ASSOCIATED FILES |
---|
[1033] | 21 | ! --> thermcell_main_mars.F90 |
---|
[1032] | 22 | ! --> thermcell_dqup.F90 |
---|
| 23 | ! --> comtherm_h.F90 |
---|
| 24 | ! ----------------------------------------------------------------------- |
---|
| 25 | ! Reference paper: |
---|
| 26 | ! A. Colaïtis, A. Spiga, F. Hourdin, C. Rio, F. Forget, and E. Millour. |
---|
| 27 | ! A thermal plume model for the Martian convective boundary layer. |
---|
| 28 | ! Journal of Geophysical Research (Planets), 118:1468-1487, July 2013. |
---|
[1033] | 29 | ! http://dx.doi.org/10.1002/jgre.20104 |
---|
| 30 | ! http://arxiv.org/abs/1306.6215 |
---|
[1032] | 31 | ! ----------------------------------------------------------------------- |
---|
| 32 | ! Reference paper for terrestrial plume model: |
---|
| 33 | ! C. Rio and F. Hourdin. |
---|
| 34 | ! A thermal plume model for the convective boundary layer : Representation of cumulus clouds. |
---|
| 35 | ! Journal of the Atmospheric Sciences, 65:407-425, 2008. |
---|
| 36 | ! ----------------------------------------------------------------------- |
---|
| 37 | |
---|
| 38 | SUBROUTINE calltherm_interface (ngrid,nlayer,nq, & |
---|
| 39 | & tracer,igcm_co2, & |
---|
[652] | 40 | & zzlev,zzlay, & |
---|
[161] | 41 | & ptimestep,pu,pv,pt,pq,pdu,pdv,pdt,pdq,q2, & |
---|
[185] | 42 | & pplay,pplev,pphi,zpopsk, & |
---|
[660] | 43 | & pdu_th,pdv_th,pdt_th,pdq_th,lmax,zmaxth,pbl_dtke, & |
---|
| 44 | & pdhdif,hfmax,wstar,sensibFlux) |
---|
[161] | 45 | |
---|
[1032] | 46 | use comtherm_h |
---|
[1036] | 47 | use tracer_mod, only: nqmx,noms |
---|
[1032] | 48 | implicit none |
---|
[161] | 49 | |
---|
[1032] | 50 | ! SHARED VARIABLES. This needs adaptations in another climate model. |
---|
[1028] | 51 | #include "comcstfi.h" !contains physical constant values such as |
---|
| 52 | ! "g" : gravitational acceleration (m.s-2) |
---|
| 53 | ! "r" : recuced gas constant (J.K-1.mol-1) |
---|
| 54 | ! "cpp" : specific heat of the atmosphere (J.kg-1.K-1) |
---|
[185] | 55 | |
---|
[161] | 56 | !-------------------------------------------------------- |
---|
[342] | 57 | ! Input Variables |
---|
[161] | 58 | !-------------------------------------------------------- |
---|
| 59 | |
---|
[1032] | 60 | INTEGER, INTENT(IN) :: ngrid ! number of horizontal grid points |
---|
| 61 | INTEGER, INTENT(IN) :: nlayer ! number of vertical grid points |
---|
| 62 | INTEGER, INTENT(IN) :: nq ! number of tracer species |
---|
[1028] | 63 | REAL, INTENT(IN) :: ptimestep !timestep (s) |
---|
[1032] | 64 | REAL, INTENT(IN) :: pplev(ngrid,nlayer+1) !intermediate pressure levels (Pa) |
---|
| 65 | REAL, INTENT(IN) :: pplay(ngrid,nlayer) !Pressure at the middle of the layers (Pa) |
---|
| 66 | REAL, INTENT(IN) :: pphi(ngrid,nlayer) !Geopotential at the middle of the layers (m2s-2) |
---|
| 67 | REAL, INTENT(IN) :: pu(ngrid,nlayer),pv(ngrid,nlayer) !u,v components of the wind (ms-1) |
---|
| 68 | REAL, INTENT(IN) :: pt(ngrid,nlayer),pq(ngrid,nlayer,nq)!temperature (K) and tracer concentration (kg/kg) |
---|
| 69 | REAL, INTENT(IN) :: zzlay(ngrid,nlayer) ! altitude at the middle of the layers |
---|
| 70 | REAL, INTENT(IN) :: zzlev(ngrid,nlayer+1) ! altitude at layer boundaries |
---|
| 71 | LOGICAL, INTENT(IN) :: tracer ! =.true. if tracers are present and to be transported |
---|
| 72 | INTEGER, INTENT(IN) :: igcm_co2 ! index of the CO2 tracer in mixing ratio array |
---|
| 73 | ! --> 0 if no tracer is CO2 (or no tracer at all) |
---|
| 74 | ! --> this prepares special treatment for polar night mixing |
---|
| 75 | ! (see thermcell_main_mars) |
---|
| 76 | REAL, INTENT(IN) :: pdu(ngrid,nlayer),pdv(ngrid,nlayer) ! wind velocity change from routines called |
---|
[1028] | 77 | ! before thermals du/dt (m/s/s) |
---|
[1032] | 78 | REAL, INTENT(IN) :: pdq(ngrid,nlayer,nq) ! tracer concentration change from routines called |
---|
[1028] | 79 | ! before thermals dq/dt (kg/kg/s) |
---|
[1032] | 80 | REAL, INTENT(IN) :: pdt(ngrid,nlayer) ! temperature change from routines called before thermals dT/dt (K/s) |
---|
| 81 | REAL, INTENT(IN) :: q2(ngrid,nlayer+1) ! turbulent kinetic energy |
---|
| 82 | REAL, INTENT(IN) :: zpopsk(ngrid,nlayer) ! ratio of pressure at middle of layer to surface pressure, |
---|
[1028] | 83 | ! to the power r/cp, i.e. zpopsk=(pplay(ig,l)/pplev(ig,1))**rcp |
---|
[1032] | 84 | REAL, INTENT(IN) :: pdhdif(ngrid,nlayer) ! potential temperature change from turbulent diffusion scheme dT/dt (K/s) |
---|
| 85 | REAL, INTENT(IN) :: sensibFlux(ngrid) ! sensible heat flux computed from surface layer scheme |
---|
[161] | 86 | |
---|
| 87 | !-------------------------------------------------------- |
---|
[342] | 88 | ! Output Variables |
---|
[161] | 89 | !-------------------------------------------------------- |
---|
| 90 | |
---|
[1032] | 91 | REAL, INTENT(OUT) :: pdu_th(ngrid,nlayer) ! wind velocity change from thermals du/dt (m/s/s) |
---|
| 92 | REAL, INTENT(OUT) :: pdv_th(ngrid,nlayer) ! wind velocity change from thermals dv/dt (m/s/s) |
---|
| 93 | REAL, INTENT(OUT) :: pdt_th(ngrid,nlayer) ! temperature change from thermals dT/dt (K/s) |
---|
| 94 | REAL, INTENT(OUT) :: pdq_th(ngrid,nlayer,nq) ! tracer change from thermals dq/dt (kg/kg/s) |
---|
| 95 | INTEGER, INTENT(OUT) :: lmax(ngrid) ! layer number reacher by thermals in grid point |
---|
| 96 | REAL, INTENT(OUT) :: zmaxth(ngrid) ! equivalent to lmax, but in (m), interpolated |
---|
| 97 | REAL, INTENT(OUT) :: pbl_dtke(ngrid,nlayer+1) ! turbulent kinetic energy change from thermals dtke/dt |
---|
| 98 | REAL, INTENT(OUT) :: wstar(ngrid) ! free convection velocity (m/s) |
---|
[161] | 99 | |
---|
| 100 | !-------------------------------------------------------- |
---|
[342] | 101 | ! Thermals local variables |
---|
[161] | 102 | !-------------------------------------------------------- |
---|
[1032] | 103 | REAL zu(ngrid,nlayer), zv(ngrid,nlayer) |
---|
| 104 | REAL zt(ngrid,nlayer) |
---|
| 105 | REAL d_t_ajs(ngrid,nlayer) |
---|
| 106 | REAL d_u_ajs(ngrid,nlayer), d_q_ajs(ngrid,nlayer,nq) |
---|
| 107 | REAL d_v_ajs(ngrid,nlayer) |
---|
| 108 | REAL fm_therm(ngrid,nlayer+1), entr_therm(ngrid,nlayer) |
---|
| 109 | REAL detr_therm(ngrid,nlayer),detrmod(ngrid,nlayer) |
---|
| 110 | REAL zw2(ngrid,nlayer+1) |
---|
| 111 | REAL fraca(ngrid,nlayer+1),zfraca(ngrid,nlayer+1) |
---|
| 112 | REAL q_therm(ngrid,nlayer), pq_therm(ngrid,nlayer,nq) |
---|
| 113 | REAL q2_therm(ngrid,nlayer), dq2_therm(ngrid,nlayer) |
---|
| 114 | REAL lmax_real(ngrid) |
---|
| 115 | REAL masse(ngrid,nlayer) |
---|
[1028] | 116 | |
---|
[660] | 117 | INTEGER l,ig,iq,ii(1),k |
---|
[628] | 118 | CHARACTER (LEN=20) modname |
---|
[161] | 119 | |
---|
[342] | 120 | !-------------------------------------------------------- |
---|
| 121 | ! Local variables for sub-timestep |
---|
| 122 | !-------------------------------------------------------- |
---|
[161] | 123 | |
---|
[1032] | 124 | REAL d_t_the(ngrid,nlayer), d_q_the(ngrid,nlayer,nq) |
---|
[561] | 125 | INTEGER isplit |
---|
[342] | 126 | REAL fact |
---|
[1032] | 127 | REAL zfm_therm(ngrid,nlayer+1),zdt |
---|
| 128 | REAL zentr_therm(ngrid,nlayer),zdetr_therm(ngrid,nlayer) |
---|
| 129 | REAL zheatFlux(ngrid,nlayer) |
---|
| 130 | REAL zheatFlux_down(ngrid,nlayer) |
---|
| 131 | REAL zbuoyancyOut(ngrid,nlayer) |
---|
| 132 | REAL zbuoyancyEst(ngrid,nlayer) |
---|
| 133 | REAL zzw2(ngrid,nlayer+1) |
---|
| 134 | REAL zmax(ngrid) |
---|
[1212] | 135 | INTEGER ndt,limz |
---|
[342] | 136 | |
---|
| 137 | !-------------------------------------------------------- |
---|
| 138 | ! Diagnostics |
---|
| 139 | !-------------------------------------------------------- |
---|
| 140 | |
---|
[1032] | 141 | REAL heatFlux(ngrid,nlayer) |
---|
| 142 | REAL heatFlux_down(ngrid,nlayer) |
---|
| 143 | REAL buoyancyOut(ngrid,nlayer) |
---|
| 144 | REAL buoyancyEst(ngrid,nlayer) |
---|
| 145 | REAL hfmax(ngrid),wmax(ngrid) |
---|
| 146 | REAL pbl_teta(ngrid),dteta(ngrid,nlayer) |
---|
| 147 | REAL rpdhd(ngrid,nlayer) |
---|
| 148 | REAL wtdif(ngrid,nlayer),rho(ngrid,nlayer) |
---|
| 149 | REAL wtth(ngrid,nlayer) |
---|
[161] | 150 | |
---|
| 151 | ! ********************************************************************** |
---|
[1028] | 152 | ! Initializations |
---|
| 153 | ! ********************************************************************** |
---|
| 154 | |
---|
[621] | 155 | lmax(:)=0 |
---|
[161] | 156 | pdu_th(:,:)=0. |
---|
| 157 | pdv_th(:,:)=0. |
---|
| 158 | pdt_th(:,:)=0. |
---|
| 159 | entr_therm(:,:)=0. |
---|
| 160 | detr_therm(:,:)=0. |
---|
| 161 | q2_therm(:,:)=0. |
---|
| 162 | dq2_therm(:,:)=0. |
---|
| 163 | pbl_dtke(:,:)=0. |
---|
| 164 | fm_therm(:,:)=0. |
---|
| 165 | zw2(:,:)=0. |
---|
| 166 | fraca(:,:)=0. |
---|
[512] | 167 | zfraca(:,:)=0. |
---|
[161] | 168 | if (tracer) then |
---|
| 169 | pdq_th(:,:,:)=0. |
---|
| 170 | end if |
---|
[342] | 171 | d_t_ajs(:,:)=0. |
---|
| 172 | d_u_ajs(:,:)=0. |
---|
| 173 | d_v_ajs(:,:)=0. |
---|
| 174 | d_q_ajs(:,:,:)=0. |
---|
| 175 | heatFlux(:,:)=0. |
---|
| 176 | heatFlux_down(:,:)=0. |
---|
| 177 | buoyancyOut(:,:)=0. |
---|
| 178 | buoyancyEst(:,:)=0. |
---|
| 179 | zmaxth(:)=0. |
---|
| 180 | lmax_real(:)=0. |
---|
[161] | 181 | |
---|
| 182 | |
---|
[342] | 183 | ! ********************************************************************** |
---|
[1028] | 184 | ! Preparing inputs for the thermals: increment tendancies |
---|
| 185 | ! from other subroutines called before the thermals model |
---|
[342] | 186 | ! ********************************************************************** |
---|
[161] | 187 | |
---|
[1028] | 188 | zu(:,:)=pu(:,:)+pdu(:,:)*ptimestep ! u-component of wind velocity |
---|
| 189 | zv(:,:)=pv(:,:)+pdv(:,:)*ptimestep ! v-component of wind velocity |
---|
| 190 | zt(:,:)=pt(:,:)+pdt(:,:)*ptimestep ! temperature |
---|
[161] | 191 | |
---|
[1028] | 192 | pq_therm(:,:,:)=0. ! tracer concentration |
---|
[161] | 193 | |
---|
[342] | 194 | if(qtransport_thermals) then |
---|
[1028] | 195 | if(tracer) then ! tracer is a logical that is true if tracer must be transported in the GCM physics |
---|
| 196 | pq_therm(:,:,:)=pq(:,:,:)+pdq(:,:,:)*ptimestep ! tracer concentration |
---|
[342] | 197 | endif |
---|
| 198 | endif |
---|
[161] | 199 | |
---|
[1028] | 200 | |
---|
[544] | 201 | IF(dtke_thermals) THEN |
---|
[1032] | 202 | DO l=1,nlayer |
---|
[544] | 203 | q2_therm(:,l)=0.5*(q2(:,l)+q2(:,l+1)) |
---|
| 204 | ENDDO |
---|
| 205 | ENDIF |
---|
[342] | 206 | |
---|
| 207 | ! ********************************************************************** |
---|
[1032] | 208 | ! --> CALLTHERM |
---|
[342] | 209 | ! SUB-TIMESTEP LOOP |
---|
| 210 | ! ********************************************************************** |
---|
| 211 | |
---|
[1028] | 212 | zdt=ptimestep/REAL(nsplit_thermals) !subtimestep |
---|
[342] | 213 | |
---|
[1028] | 214 | DO isplit=1,nsplit_thermals !temporal loop on the subtimestep |
---|
[342] | 215 | |
---|
| 216 | ! Initialization of intermediary variables |
---|
| 217 | |
---|
| 218 | zzw2(:,:)=0. |
---|
| 219 | zmax(:)=0. |
---|
[621] | 220 | lmax(:)=0 |
---|
[628] | 221 | |
---|
[1032] | 222 | if (nq .ne. 0 .and. igcm_co2 .ne. 0) then !initialize co2 tracer tendancy |
---|
| 223 | d_q_the(:,:,igcm_co2)=0. |
---|
[161] | 224 | endif |
---|
| 225 | |
---|
[1028] | 226 | ! CALL to main thermal routine |
---|
[1032] | 227 | CALL thermcell_main_mars(ngrid,nlayer,nq & |
---|
| 228 | & ,tracer,igcm_co2 & |
---|
| 229 | & ,zdt & |
---|
[342] | 230 | & ,pplay,pplev,pphi,zzlev,zzlay & |
---|
| 231 | & ,zu,zv,zt,pq_therm,q2_therm & |
---|
[1028] | 232 | & ,d_t_the,d_q_the & |
---|
[1212] | 233 | & ,zfm_therm,zentr_therm,zdetr_therm,lmax,zmax,limz & |
---|
[342] | 234 | & ,zzw2,fraca,zpopsk & |
---|
[1028] | 235 | & ,zheatFlux,zheatFlux_down & |
---|
[342] | 236 | & ,zbuoyancyOut,zbuoyancyEst) |
---|
[161] | 237 | |
---|
[342] | 238 | fact=1./REAL(nsplit_thermals) |
---|
[161] | 239 | |
---|
[1028] | 240 | ! Update thermals tendancies |
---|
| 241 | d_t_the(:,:)=d_t_the(:,:)*ptimestep*fact !temperature |
---|
[1032] | 242 | if (igcm_co2 .ne. 0) then |
---|
| 243 | d_q_the(:,:,igcm_co2)=d_q_the(:,:,igcm_co2)*ptimestep*fact !co2 mass mixing ratio |
---|
[508] | 244 | endif |
---|
[1028] | 245 | zmaxth(:)=zmaxth(:)+zmax(:)*fact !thermals height |
---|
| 246 | lmax_real(:)=lmax_real(:)+float(lmax(:))*fact !thermals height index |
---|
| 247 | fm_therm(:,:)=fm_therm(:,:) & !upward mass flux |
---|
[342] | 248 | & +zfm_therm(:,:)*fact |
---|
[1028] | 249 | entr_therm(:,:)=entr_therm(:,:) & !entrainment mass flux |
---|
[342] | 250 | & +zentr_therm(:,:)*fact |
---|
[1028] | 251 | detr_therm(:,:)=detr_therm(:,:) & !detrainment mass flux |
---|
[342] | 252 | & +zdetr_therm(:,:)*fact |
---|
[1028] | 253 | zfraca(:,:)=zfraca(:,:) + fraca(:,:)*fact !updraft fractional coverage |
---|
| 254 | heatFlux(:,:)=heatFlux(:,:) & !upward heat flux |
---|
[342] | 255 | & +zheatFlux(:,:)*fact |
---|
[1028] | 256 | heatFlux_down(:,:)=heatFlux_down(:,:) & !downward heat flux |
---|
[342] | 257 | & +zheatFlux_down(:,:)*fact |
---|
[1028] | 258 | buoyancyOut(:,:)=buoyancyOut(:,:) & !plume final buoyancy |
---|
[508] | 259 | & +zbuoyancyOut(:,:)*fact |
---|
[1028] | 260 | buoyancyEst(:,:)=buoyancyEst(:,:) & !plume estimated buoyancy used for vertical velocity computation |
---|
[508] | 261 | & +zbuoyancyEst(:,:)*fact |
---|
[1028] | 262 | zw2(:,:)=zw2(:,:) + zzw2(:,:)*fact !vertical velocity |
---|
[342] | 263 | |
---|
[1028] | 264 | ! Save tendancies |
---|
| 265 | d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:) !temperature tendancy (delta T) |
---|
[1032] | 266 | if (igcm_co2 .ne. 0) then |
---|
| 267 | d_q_ajs(:,:,igcm_co2)=d_q_ajs(:,:,igcm_co2)+d_q_the(:,:,igcm_co2) !tracer tendancy (delta q) |
---|
[508] | 268 | endif |
---|
[342] | 269 | |
---|
[1028] | 270 | ! Increment temperature and co2 concentration for next pass in subtimestep loop |
---|
| 271 | zt(:,:) = zt(:,:) + d_t_the(:,:) !temperature |
---|
[1032] | 272 | if (igcm_co2 .ne. 0) then |
---|
| 273 | pq_therm(:,:,igcm_co2) = & |
---|
| 274 | & pq_therm(:,:,igcm_co2) + d_q_the(:,:,igcm_co2) !co2 tracer |
---|
[508] | 275 | endif |
---|
[342] | 276 | |
---|
| 277 | |
---|
| 278 | ENDDO ! isplit |
---|
| 279 | !**************************************************************** |
---|
| 280 | |
---|
| 281 | ! Now that we have computed total entrainment and detrainment, we can |
---|
[1028] | 282 | ! advect u, v, and q in thermals. (potential temperature and co2 MMR |
---|
| 283 | ! have already been advected in thermcell_main because they are coupled |
---|
| 284 | ! to the determination of the thermals caracteristics). This is done |
---|
| 285 | ! separatly because u,v, and q are not used in thermcell_main for |
---|
[342] | 286 | ! any thermals-related computation : they are purely passive. |
---|
| 287 | |
---|
| 288 | ! mass of cells |
---|
[1032] | 289 | do l=1,nlayer |
---|
[342] | 290 | masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g |
---|
| 291 | enddo |
---|
| 292 | |
---|
[1028] | 293 | ! recompute detrainment mass flux from entrainment and upward mass flux |
---|
| 294 | ! this ensure mass flux conservation |
---|
[628] | 295 | detrmod(:,:)=0. |
---|
[1212] | 296 | do l=1,limz |
---|
[1032] | 297 | do ig=1,ngrid |
---|
[628] | 298 | detrmod(ig,l)=fm_therm(ig,l)-fm_therm(ig,l+1) & |
---|
| 299 | & +entr_therm(ig,l) |
---|
| 300 | if (detrmod(ig,l).lt.0.) then |
---|
| 301 | entr_therm(ig,l)=entr_therm(ig,l)-detrmod(ig,l) |
---|
| 302 | detrmod(ig,l)=0. |
---|
| 303 | endif |
---|
| 304 | enddo |
---|
| 305 | enddo |
---|
[1028] | 306 | |
---|
| 307 | ! u component of wind velocity advection in thermals |
---|
| 308 | ! result is a derivative (d_u_ajs in m/s/s) |
---|
[628] | 309 | ndt=10 |
---|
[1032] | 310 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
[628] | 311 | & ,fm_therm,entr_therm,detrmod, & |
---|
[1212] | 312 | & masse,zu,d_u_ajs,ndt,limz) |
---|
[342] | 313 | |
---|
[1028] | 314 | ! v component of wind velocity advection in thermals |
---|
| 315 | ! result is a derivative (d_v_ajs in m/s/s) |
---|
[1032] | 316 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
[628] | 317 | & ,fm_therm,entr_therm,detrmod, & |
---|
[1212] | 318 | & masse,zv,d_v_ajs,ndt,limz) |
---|
[628] | 319 | |
---|
[1028] | 320 | ! non co2 tracers advection in thermals |
---|
| 321 | ! result is a derivative (d_q_ajs in kg/kg/s) |
---|
| 322 | |
---|
[1032] | 323 | if (nq .ne. 0.) then |
---|
| 324 | DO iq=1,nq |
---|
| 325 | if (iq .ne. igcm_co2) then |
---|
| 326 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
[628] | 327 | & ,fm_therm,entr_therm,detrmod, & |
---|
[1212] | 328 | & masse,pq_therm(:,:,iq),d_q_ajs(:,:,iq),ndt,limz) |
---|
[508] | 329 | endif |
---|
[342] | 330 | ENDDO |
---|
| 331 | endif |
---|
| 332 | |
---|
[1028] | 333 | ! tke advection in thermals |
---|
| 334 | ! result is a tendancy (d_u_ajs in J) |
---|
[544] | 335 | if (dtke_thermals) then |
---|
[1032] | 336 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
[628] | 337 | & ,fm_therm,entr_therm,detrmod, & |
---|
[1212] | 338 | & masse,q2_therm,dq2_therm,ndt,limz) |
---|
[544] | 339 | endif |
---|
| 340 | |
---|
[1028] | 341 | ! compute wmax for diagnostics |
---|
[1032] | 342 | DO ig=1,ngrid |
---|
[342] | 343 | wmax(ig)=MAXVAL(zw2(ig,:)) |
---|
| 344 | ENDDO |
---|
| 345 | |
---|
| 346 | ! ********************************************************************** |
---|
| 347 | ! ********************************************************************** |
---|
| 348 | ! ********************************************************************** |
---|
| 349 | ! CALLTHERM END |
---|
| 350 | ! ********************************************************************** |
---|
| 351 | ! ********************************************************************** |
---|
| 352 | ! ********************************************************************** |
---|
| 353 | |
---|
| 354 | |
---|
| 355 | ! ********************************************************************** |
---|
| 356 | ! Preparing outputs |
---|
| 357 | ! ********************************************************************** |
---|
| 358 | |
---|
[1212] | 359 | do l=1,limz |
---|
[628] | 360 | pdu_th(:,l)=d_u_ajs(:,l) |
---|
| 361 | pdv_th(:,l)=d_v_ajs(:,l) |
---|
| 362 | enddo |
---|
[342] | 363 | |
---|
[1028] | 364 | ! if tracers are transported in thermals, update output variables, else these are 0. |
---|
[161] | 365 | if(qtransport_thermals) then |
---|
[342] | 366 | if(tracer) then |
---|
[1032] | 367 | do iq=1,nq |
---|
| 368 | if (iq .ne. igcm_co2) then |
---|
[1212] | 369 | do l=1,limz |
---|
[1028] | 370 | pdq_th(:,l,iq)=d_q_ajs(:,l,iq) !non-co2 tracers d_q_ajs are dq/dt (kg/kg/s) |
---|
[628] | 371 | enddo |
---|
[625] | 372 | else |
---|
[1212] | 373 | do l=1,limz |
---|
[1028] | 374 | pdq_th(:,l,iq)=d_q_ajs(:,l,iq)/ptimestep !co2 tracer d_q_ajs is dq (kg/kg) |
---|
[628] | 375 | enddo |
---|
[625] | 376 | endif |
---|
| 377 | enddo |
---|
[342] | 378 | endif |
---|
[161] | 379 | endif |
---|
| 380 | |
---|
[1028] | 381 | ! if tke is transported in thermals, update output variable, else this is 0. |
---|
[544] | 382 | IF(dtke_thermals) THEN |
---|
[1032] | 383 | DO l=2,nlayer |
---|
[544] | 384 | pbl_dtke(:,l)=0.5*(dq2_therm(:,l-1)+dq2_therm(:,l)) |
---|
| 385 | ENDDO |
---|
| 386 | |
---|
| 387 | pbl_dtke(:,1)=0.5*dq2_therm(:,1) |
---|
[1032] | 388 | pbl_dtke(:,nlayer+1)=0. |
---|
[544] | 389 | ENDIF |
---|
[161] | 390 | |
---|
[1028] | 391 | ! update output variable for temperature. d_t_ajs is delta T in (K), pdt_th is dT/dt in (K/s) |
---|
[1212] | 392 | do l=1,limz |
---|
[628] | 393 | pdt_th(:,l)=d_t_ajs(:,l)/ptimestep |
---|
| 394 | enddo |
---|
[342] | 395 | |
---|
[499] | 396 | |
---|
[342] | 397 | ! ********************************************************************** |
---|
[1028] | 398 | ! SURFACE LAYER INTERFACE |
---|
| 399 | ! Compute the free convection velocity w* scale for surface layer gustiness |
---|
| 400 | ! speed parameterization. The computed value of w* will be used at the next |
---|
| 401 | ! timestep to modify surface-atmosphere exchange fluxes, because the surface |
---|
| 402 | ! layer scheme and diffusion are called BEFORE the thermals. (outside of these |
---|
| 403 | ! routines) |
---|
[499] | 404 | ! ********************************************************************** |
---|
| 405 | |
---|
| 406 | ! Potential temperature gradient |
---|
| 407 | |
---|
[1032] | 408 | dteta(:,nlayer)=0. |
---|
| 409 | DO l=1,nlayer-1 |
---|
| 410 | DO ig=1, ngrid |
---|
[499] | 411 | dteta(ig,l) = ((zt(ig,l+1)-zt(ig,l))/zpopsk(ig,l)) & |
---|
| 412 | & /(zzlay(ig,l+1)-zzlay(ig,l)) |
---|
| 413 | ENDDO |
---|
| 414 | ENDDO |
---|
| 415 | |
---|
[1028] | 416 | ! Computation of the PBL mixed layer temperature |
---|
[499] | 417 | |
---|
[1032] | 418 | DO ig=1, ngrid |
---|
[499] | 419 | ii=MINLOC(abs(dteta(ig,1:lmax(ig)))) |
---|
| 420 | pbl_teta(ig) = zt(ig,ii(1))/zpopsk(ig,ii(1)) |
---|
| 421 | ENDDO |
---|
| 422 | |
---|
[1028] | 423 | ! In order to have an accurate w*, we must add the heat flux from the |
---|
| 424 | ! diffusion scheme to the computation of the maximum heat flux hfmax |
---|
| 425 | ! Here pdhdif is the potential temperature change from the diffusion |
---|
| 426 | ! scheme (Mellor and Yamada, see paper section 6, paragraph 57) |
---|
[660] | 427 | |
---|
| 428 | ! compute rho as it is after the diffusion |
---|
| 429 | |
---|
| 430 | rho(:,:)=pplay(:,:) & |
---|
| 431 | & /(r*(pt(:,:)+pdhdif(:,:)*zpopsk(:,:)*ptimestep)) |
---|
| 432 | |
---|
| 433 | ! integrate -rho*pdhdif |
---|
| 434 | |
---|
| 435 | rpdhd(:,:)=0. |
---|
| 436 | |
---|
[1032] | 437 | DO ig=1,ngrid |
---|
[660] | 438 | DO l=1,lmax(ig) |
---|
| 439 | rpdhd(ig,l)=0. |
---|
| 440 | DO k=1,l |
---|
| 441 | rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*pdhdif(ig,k)* & |
---|
| 442 | & (zzlev(ig,k+1)-zzlev(ig,k)) |
---|
| 443 | ENDDO |
---|
| 444 | rpdhd(ig,l)=rpdhd(ig,l)-sensibFlux(ig)/cpp |
---|
| 445 | ENDDO |
---|
| 446 | ENDDO |
---|
| 447 | |
---|
[1028] | 448 | ! compute w'theta' (vertical turbulent flux of temperature) from |
---|
| 449 | ! the diffusion scheme |
---|
[660] | 450 | |
---|
| 451 | wtdif(:,:)=rpdhd(:,:)/rho(:,:) |
---|
| 452 | |
---|
[1028] | 453 | ! Now we compute the contribution of the thermals to w'theta': |
---|
[660] | 454 | ! compute rho as it is after the thermals |
---|
| 455 | |
---|
| 456 | rho(:,:)=pplay(:,:)/(r*(zt(:,:))) |
---|
[1028] | 457 | |
---|
[660] | 458 | ! integrate -rho*pdhdif |
---|
| 459 | |
---|
[1032] | 460 | DO ig=1,ngrid |
---|
[660] | 461 | DO l=1,lmax(ig) |
---|
| 462 | rpdhd(ig,l)=0. |
---|
| 463 | DO k=1,l |
---|
| 464 | rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*(pdt_th(ig,k)/zpopsk(ig,k))* & |
---|
| 465 | & (zzlev(ig,k+1)-zzlev(ig,k)) |
---|
| 466 | ENDDO |
---|
| 467 | rpdhd(ig,l)=rpdhd(ig,l)+ & |
---|
| 468 | & rho(ig,1)*(heatFlux(ig,1)+heatFlux_down(ig,1)) |
---|
| 469 | ENDDO |
---|
| 470 | ENDDO |
---|
[1032] | 471 | rpdhd(:,nlayer)=0. |
---|
[660] | 472 | |
---|
| 473 | ! compute w'teta' from thermals |
---|
| 474 | |
---|
| 475 | wtth(:,:)=rpdhd(:,:)/rho(:,:) |
---|
| 476 | |
---|
[1028] | 477 | ! Add vertical turbulent heat fluxes from the thermals and the diffusion scheme |
---|
| 478 | ! and compute the maximum |
---|
[660] | 479 | |
---|
[1032] | 480 | DO ig=1,ngrid |
---|
[660] | 481 | hfmax(ig)=MAXVAL(wtth(ig,:)+wtdif(ig,:)) |
---|
| 482 | ENDDO |
---|
[1028] | 483 | |
---|
| 484 | ! Finally we can compute the free convection velocity scale |
---|
[499] | 485 | ! We follow Spiga et. al 2010 (QJRMS) |
---|
| 486 | ! ------------ |
---|
| 487 | |
---|
[1032] | 488 | DO ig=1, ngrid |
---|
[499] | 489 | IF (zmax(ig) .gt. 0.) THEN |
---|
| 490 | wstar(ig)=(g*zmaxth(ig)*hfmax(ig)/pbl_teta(ig))**(1./3.) |
---|
| 491 | ELSE |
---|
| 492 | wstar(ig)=0. |
---|
| 493 | ENDIF |
---|
| 494 | ENDDO |
---|
| 495 | |
---|
[161] | 496 | END |
---|