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