[3175] | 1 | SUBROUTINE ch4surf(ngrid,nlayer,nq,ptimestep,tsurf,pdtsurf, |
---|
| 2 | & pplev,pdpsurf,pq,pdq,pqsurf,pdqsurf,pdqch4,pdqsch4) |
---|
| 3 | |
---|
| 4 | use planet_h |
---|
| 5 | use comgeomfi_h |
---|
| 6 | IMPLICIT NONE |
---|
| 7 | |
---|
| 8 | c---------------- |
---|
| 9 | c declarations: |
---|
| 10 | c ------------- |
---|
| 11 | |
---|
| 12 | #include "dimensions.h" |
---|
| 13 | #include "dimphys.h" |
---|
| 14 | #include "comcstfi.h" |
---|
| 15 | #include "comsaison.h" |
---|
| 16 | #include "callkeys.h" |
---|
| 17 | #include "surfdat.h" |
---|
| 18 | #include "tracer.h" |
---|
| 19 | !#include "planete.h" |
---|
| 20 | |
---|
| 21 | ! Routine for nogcm : sublimation/condensation scheme at the surface |
---|
| 22 | ! Output : tendancy for methane mixing ratio and surface reservoir : |
---|
| 23 | ! pdqch4, pdqs_ch4 |
---|
| 24 | |
---|
| 25 | !----------------------------------------------------------------------- |
---|
| 26 | ! Arguments |
---|
| 27 | |
---|
| 28 | |
---|
| 29 | INTEGER ngrid,nlayer,nq |
---|
| 30 | REAL ptimestep |
---|
| 31 | INTEGER ig,iq |
---|
| 32 | |
---|
| 33 | ! input : |
---|
| 34 | REAL tsurf(ngrid) |
---|
| 35 | REAL pplev(ngrid,nlayer+1) |
---|
| 36 | REAL pdpsurf(ngrid) |
---|
| 37 | REAL pq(ngrid,nlayer,nq) |
---|
| 38 | REAL pdq(ngrid,nlayer,nq) |
---|
| 39 | REAL pqsurf(ngrid,nq) |
---|
| 40 | REAL pdqsurf(ngrid,nq) |
---|
| 41 | REAL pdtsurf(ngrid) |
---|
| 42 | |
---|
| 43 | ! Output |
---|
| 44 | REAL pdqch4(ngrid) |
---|
| 45 | REAL pdqsch4(ngrid) |
---|
| 46 | |
---|
| 47 | ! local |
---|
| 48 | REAL qsat(ngridmx) |
---|
| 49 | REAL zpsrf(ngridmx) |
---|
| 50 | REAL zq_ch4(ngridmx) |
---|
| 51 | REAL zq_n2surf(ngridmx) |
---|
| 52 | REAL ztsurf(ngridmx) |
---|
| 53 | REAL gamm(ngridmx) ! activity coefficient |
---|
| 54 | REAL rho,u,v,uv,z00,cdrag,alt |
---|
| 55 | REAL vonk ! Von Karman Constant |
---|
| 56 | SAVE vonk |
---|
| 57 | DATA vonk/0.4/ |
---|
| 58 | |
---|
| 59 | ! Calculation of turbulent flux : F=rho*cdrag*uv*(qsat-zq) |
---|
| 60 | |
---|
| 61 | ! Calcul de cdrag |
---|
| 62 | alt=5. ! m |
---|
| 63 | !z00=1.e-2 ! rugosity |
---|
| 64 | z00=z0 |
---|
| 65 | cdrag=(vonk/log(alt/z00))**2 |
---|
| 66 | |
---|
| 67 | u=6. ! 6 |
---|
| 68 | v=3. ! 3 |
---|
| 69 | uv=sqrt(u**2+v**2) |
---|
| 70 | |
---|
| 71 | pdqsch4(:)=0. |
---|
| 72 | pdqch4(:)=0. |
---|
| 73 | |
---|
| 74 | !! Update CH4, pressure |
---|
| 75 | DO ig=1,ngridmx |
---|
| 76 | zpsrf(ig)=pplev(ig,1) |
---|
| 77 | zq_ch4(ig)=pq(ig,1,igcm_ch4_gas) |
---|
| 78 | & + pdq(ig,1,igcm_ch4_gas)*ptimestep |
---|
| 79 | zq_n2surf(ig)=pqsurf(ig,igcm_n2) |
---|
| 80 | & + pdqsurf(ig,igcm_n2)*ptimestep |
---|
| 81 | ztsurf(ig)=tsurf(ig) !+pdtsurf(ig)*ptimestep |
---|
| 82 | ENDDO |
---|
| 83 | |
---|
| 84 | !! Get qsat for CH4 |
---|
| 85 | call methanesat(ngridmx,ztsurf,zpsrf,qsat,zq_n2surf(:)) |
---|
| 86 | |
---|
| 87 | !! Dayfrac: Fraction of the daytime where we do not condense CH4 in N2 |
---|
| 88 | ! corresponds to cold layer of N2 pushing CH4 or depleting CH4 near the surface |
---|
| 89 | ! By default, we do not take this into account : dayfrac=0. We condense all the time |
---|
| 90 | ! dayfrac = 1 : we do not condense CH4 in N2 during daytime |
---|
| 91 | ! dayfrac = 2 : all saturated, so we never condense CH4 in N2 |
---|
| 92 | ! dayfrac = 0.5 : we do not condense CH4 during half of the daytime (we condense during all night + half of daytime) |
---|
| 93 | |
---|
| 94 | !! Loop |
---|
| 95 | DO ig=1,ngridmx |
---|
| 96 | |
---|
| 97 | !! Take into account activity coefficient |
---|
| 98 | gamm(ig)=499.9-21.8*ztsurf(ig)+0.249*ztsurf(ig)**2 |
---|
| 99 | & -1.3*(zq_ch4(ig)* |
---|
| 100 | & mmol(igcm_n2)/mmol(igcm_ch4_gas)*100.-0.6)/0.3 |
---|
| 101 | gamm(ig)=max(gamm(ig),1.) |
---|
| 102 | qsat(ig)=qsat(ig)*gamm(ig) |
---|
| 103 | |
---|
| 104 | rho = zpsrf(ig) / (r * tsurf(ig) ) |
---|
| 105 | !! Condensation Flux |
---|
| 106 | pdqsch4(ig)=(-rho*uv*cdrag*(qsat(ig)-zq_ch4(ig))) |
---|
| 107 | |
---|
| 108 | if (dayfrac.gt.0.and.zq_n2surf(ig).gt.thresh_non2) then |
---|
| 109 | if (dayfrac.gt.1.) then |
---|
| 110 | pdqsch4(ig)=min(pdqsch4(ig),0.) |
---|
| 111 | else |
---|
| 112 | if (pdqsch4(ig).gt.0.) then ! condensation |
---|
| 113 | pdqsch4(ig)=pdqsch4(ig)*(1.-fract(ig)*dayfrac) |
---|
| 114 | endif |
---|
| 115 | endif |
---|
| 116 | endif |
---|
| 117 | |
---|
| 118 | !! Conserve mass if reservoir depleted |
---|
| 119 | if ((-pdqsch4(ig)*ptimestep).gt. |
---|
| 120 | & (pqsurf(ig,igcm_ch4_ice))) then |
---|
| 121 | pdqsch4(ig)=-pqsurf(ig,igcm_ch4_ice)/ptimestep |
---|
| 122 | endif |
---|
| 123 | |
---|
| 124 | if (pdqsch4(ig)*ptimestep.gt.zq_ch4(ig)*zpsrf(ig)/g) then |
---|
| 125 | pdqsch4(ig)=zq_ch4(ig)/ptimestep*zpsrf(ig)/g |
---|
| 126 | endif |
---|
| 127 | |
---|
| 128 | !! Security to avoid large changes in temperatures due to |
---|
| 129 | !latent heat |
---|
| 130 | if (pdqsch4(ig)*ptimestep.gt.0.25) then |
---|
| 131 | pdqsch4(ig)=0.25/ptimestep |
---|
| 132 | endif |
---|
| 133 | if (pdqsch4(ig)*ptimestep.lt.-0.25) then |
---|
| 134 | pdqsch4(ig)=-0.25/ptimestep |
---|
| 135 | endif |
---|
| 136 | |
---|
| 137 | !! Atm tendency |
---|
| 138 | pdqch4(ig)=-pdqsch4(ig)*g/zpsrf(ig) |
---|
| 139 | |
---|
| 140 | ENDDO |
---|
| 141 | |
---|
| 142 | RETURN |
---|
| 143 | END |
---|