[3739] | 1 | MODULE co2snow_mod |
---|
| 2 | |
---|
| 3 | IMPLICIT NONE |
---|
| 4 | |
---|
| 5 | CONTAINS |
---|
[2362] | 6 | c======================================================================= |
---|
| 7 | c Program for simulate the impact of the CO2 snow fall on |
---|
| 8 | c the surface infrared emission (emissivity) |
---|
| 9 | c----------------------------------------------------------------------- |
---|
| 10 | c Algorithme: |
---|
| 11 | c 1 - Initialization |
---|
| 12 | c 2 - Compute the surface emissivity |
---|
| 13 | c----------------------------------------------------------------------- |
---|
| 14 | c Author: F. Forget |
---|
| 15 | c----------------------------------------------------------------------- |
---|
| 16 | c Reference paper: |
---|
| 17 | c Francois Forget and James B. Pollack |
---|
| 18 | c |
---|
| 19 | c "Thermal infrared observations of the condensing Martian polar |
---|
| 20 | c caps: CO2 ice temperatures and radiative budget" |
---|
| 21 | c |
---|
| 22 | c JOURNAL OF GEOPHYSICAL RESEARCH, VOL. 101, NO. E7, PAGES 16, |
---|
| 23 | c 865-16,879, JULY 25, 1996 |
---|
| 24 | c======================================================================= |
---|
[38] | 25 | |
---|
[2362] | 26 | SUBROUTINE co2snow(ngrid, nlayer, ptimestep, emisref, condsub, |
---|
| 27 | & pplev, pcondicea, pcondices, pfallice, |
---|
| 28 | & pemisurf) |
---|
| 29 | |
---|
[1047] | 30 | use surfdat_h, only: iceradius, dtemisice |
---|
[1543] | 31 | use geometry_mod, only: latitude ! grid point latitudes (rad) |
---|
[1524] | 32 | use time_phylmdz_mod, only: daysec |
---|
[2362] | 33 | |
---|
[1047] | 34 | IMPLICIT NONE |
---|
[38] | 35 | |
---|
[2362] | 36 | c======================================================================= |
---|
| 37 | c Variables |
---|
| 38 | c======================================================================= |
---|
| 39 | c Inputs: |
---|
| 40 | c ------- |
---|
| 41 | integer, intent(in) :: |
---|
| 42 | & ngrid, ! number of atmospheric columns |
---|
| 43 | & nlayer ! number of atmospheric layers |
---|
| 44 | |
---|
| 45 | real, intent(in) :: |
---|
| 46 | & ptimestep, ! timestep of the physics (s) |
---|
| 47 | & emisref(ngrid), ! grd or ice emissivity without snow |
---|
| 48 | & pplev(ngrid,nlayer+1), ! interlayer pressure (Pa) |
---|
| 49 | & pcondicea(ngrid,nlayer), ! CO2 condensation rate in layers |
---|
| 50 | ! (kg/m2/s) |
---|
| 51 | & pcondices(ngrid), ! CO2 condensation rate on the surface |
---|
| 52 | ! (kg/m2/s) |
---|
| 53 | & pfallice(ngrid,nlayer+1) ! falling CO2 ice (kg/m2/s) |
---|
[2447] | 54 | |
---|
[2362] | 55 | logical, intent(in) :: |
---|
| 56 | & condsub(ngrid) ! true if there is CO2 condensation or |
---|
| 57 | ! sublimation in the column |
---|
[38] | 58 | |
---|
[2362] | 59 | c Output: |
---|
| 60 | c ------- |
---|
| 61 | real, intent(out) :: |
---|
| 62 | & pemisurf(ngrid) ! surface emissivity |
---|
[38] | 63 | |
---|
[2362] | 64 | c local: |
---|
| 65 | c ------ |
---|
| 66 | integer :: |
---|
| 67 | & l, |
---|
| 68 | & ig, |
---|
| 69 | & icap |
---|
[38] | 70 | |
---|
[2362] | 71 | real :: |
---|
| 72 | & sumdaer, |
---|
| 73 | & zdemisurf |
---|
[38] | 74 | |
---|
[2362] | 75 | real, parameter :: |
---|
| 76 | & alpha = 0.45 |
---|
[38] | 77 | |
---|
[2362] | 78 | c saved: |
---|
| 79 | c ------ |
---|
| 80 | real, save :: |
---|
| 81 | & Kscat(2) ! Kscat: coefficient for decreasing the surface |
---|
| 82 | c ! emissivity |
---|
| 83 | c ! |
---|
| 84 | c ! Kscat = (0.001/3.) * alpha / iceradius |
---|
| 85 | c ! with 0.3 < alpha < 0.6 |
---|
| 86 | c ! |
---|
| 87 | c ! alpha set to 0.45 (coeff from emis = f (tau)) and |
---|
| 88 | c ! iceradius the mean radius of the scaterring |
---|
| 89 | c ! particles (200.e-6 < iceradius < 10.e-6) |
---|
| 90 | c ! |
---|
| 91 | c ! (2) = N and S hemispheres |
---|
[38] | 92 | |
---|
[2362] | 93 | logical, save :: |
---|
| 94 | & firstcall = .true. |
---|
[2616] | 95 | |
---|
| 96 | !$OMP THREADPRIVATE(Kscat,firstcall) |
---|
| 97 | |
---|
[2362] | 98 | c======================================================================= |
---|
| 99 | c BEGIN |
---|
| 100 | c======================================================================= |
---|
| 101 | c 1 - Initialization |
---|
| 102 | c======================================================================= |
---|
[38] | 103 | if (firstcall) then |
---|
[2362] | 104 | Kscat(1) = (0.001/3.) * alpha / iceradius(1) |
---|
| 105 | Kscat(2) = (0.001/3.) * alpha / iceradius(2) |
---|
| 106 | |
---|
| 107 | firstcall = .false. |
---|
[38] | 108 | end if |
---|
[2362] | 109 | c======================================================================= |
---|
| 110 | c 2 - Compute the surface emissivity |
---|
| 111 | c======================================================================= |
---|
| 112 | do ig = 1, ngrid |
---|
| 113 | if (condsub(ig)) then |
---|
[2447] | 114 | |
---|
[2362] | 115 | if (latitude(ig).lt.0.) then |
---|
| 116 | icap = 2 ! Southern hemisphere |
---|
| 117 | else |
---|
| 118 | icap = 1 ! Northern Hemisphere |
---|
| 119 | end if |
---|
[38] | 120 | |
---|
[2362] | 121 | ! compute zdemisurf using an integrated form for numerical |
---|
| 122 | ! safety instead |
---|
| 123 | zdemisurf = |
---|
| 124 | & (emisref(ig) - pemisurf(ig)) / (dtemisice(icap) * daysec) |
---|
| 125 | & + |
---|
| 126 | & ( emisref(ig) * |
---|
| 127 | & ( (pemisurf(ig)/emisref(ig))**(-3) + |
---|
| 128 | & 3.*Kscat(icap)*pfallice(ig,1)*ptimestep )**(-1/3.) |
---|
| 129 | & - pemisurf(ig) |
---|
| 130 | & ) / ptimestep |
---|
| 131 | pemisurf(ig) = pemisurf(ig) + zdemisurf * ptimestep |
---|
| 132 | if (pemisurf(ig).lt.0.1) then |
---|
| 133 | write(*,*)'ds co2snow: emis < 0.1 !!!' |
---|
| 134 | write(*,*)'ig =' , ig |
---|
| 135 | write(*,*)'pemisurf(ig)', pemisurf(ig) |
---|
| 136 | write(*,*)'zdemisurf*ptimestep', zdemisurf*ptimestep |
---|
| 137 | end if |
---|
| 138 | else ! if condsub(ig) is false |
---|
[38] | 139 | pemisurf(ig) = emisref(ig) |
---|
[2362] | 140 | end if |
---|
| 141 | end do ! ngrid |
---|
| 142 | c======================================================================= |
---|
| 143 | c END |
---|
| 144 | c======================================================================= |
---|
[3739] | 145 | end subroutine co2snow |
---|
| 146 | |
---|
| 147 | END MODULE co2snow_mod |
---|