1 | module jthermcalc_util |
---|
2 | |
---|
3 | implicit none |
---|
4 | |
---|
5 | contains |
---|
6 | |
---|
7 | c********************************************************************** |
---|
8 | c********************************************************************** |
---|
9 | |
---|
10 | subroutine column(ig,nlayer,chemthermod,rm,nesptherm,tx,iz,zenit, |
---|
11 | $ co2colx,o2colx,o3pcolx,h2colx,h2ocolx,h2o2colx,o3colx, |
---|
12 | $ n2colx,ncolx,nocolx,cocolx,hcolx,no2colx) |
---|
13 | |
---|
14 | c nov 2002 fgg first version |
---|
15 | |
---|
16 | c********************************************************************** |
---|
17 | |
---|
18 | use tracer_mod, only: igcm_o, igcm_co2, igcm_o2, igcm_h2, |
---|
19 | & igcm_h2o_vap, igcm_h2o2, igcm_co, igcm_h, |
---|
20 | & igcm_o3, igcm_n2, igcm_n, igcm_no, igcm_no2, |
---|
21 | & mmol |
---|
22 | use param_v4_h, only: radio,gg,masa,kboltzman,n_avog |
---|
23 | |
---|
24 | implicit none |
---|
25 | |
---|
26 | |
---|
27 | c common variables and constants |
---|
28 | include 'callkeys.h' |
---|
29 | |
---|
30 | |
---|
31 | |
---|
32 | c local parameters and variables |
---|
33 | |
---|
34 | |
---|
35 | |
---|
36 | c input and output variables |
---|
37 | |
---|
38 | integer ig,nlayer |
---|
39 | integer chemthermod |
---|
40 | integer nesptherm !# of species undergoing chemistry, input |
---|
41 | real rm(nlayer,nesptherm) !densities (cm-3), input |
---|
42 | real tx(nlayer) !temperature profile, input |
---|
43 | real iz(nlayer+1) !height profile, input |
---|
44 | real zenit !SZA, input |
---|
45 | real co2colx(nlayer) !column density of CO2 (cm^-2), output |
---|
46 | real o2colx(nlayer) !column density of O2(cm^-2), output |
---|
47 | real o3pcolx(nlayer) !column density of O(3P)(cm^-2), output |
---|
48 | real h2colx(nlayer) !H2 column density (cm-2), output |
---|
49 | real h2ocolx(nlayer) !H2O column density (cm-2), output |
---|
50 | real h2o2colx(nlayer) !column density of H2O2(cm^-2), output |
---|
51 | real o3colx(nlayer) !O3 column density (cm-2), output |
---|
52 | real n2colx(nlayer) !N2 column density (cm-2), output |
---|
53 | real ncolx(nlayer) !N column density (cm-2), output |
---|
54 | real nocolx(nlayer) !NO column density (cm-2), output |
---|
55 | real cocolx(nlayer) !CO column density (cm-2), output |
---|
56 | real hcolx(nlayer) !H column density (cm-2), output |
---|
57 | real no2colx(nlayer) !NO2 column density (cm-2), output |
---|
58 | |
---|
59 | |
---|
60 | c local variables |
---|
61 | |
---|
62 | real xx |
---|
63 | real grav(nlayer) |
---|
64 | real Hco2,Ho3p,Ho2,Hh2,Hh2o,Hh2o2 |
---|
65 | real Ho3,Hn2,Hn,Hno,Hco,Hh,Hno2 |
---|
66 | |
---|
67 | real co2x(nlayer) |
---|
68 | real o2x(nlayer) |
---|
69 | real o3px(nlayer) |
---|
70 | real cox(nlayer) |
---|
71 | real hx(nlayer) |
---|
72 | real h2x(nlayer) |
---|
73 | real h2ox(nlayer) |
---|
74 | real h2o2x(nlayer) |
---|
75 | real o3x(nlayer) |
---|
76 | real n2x(nlayer) |
---|
77 | real nx(nlayer) |
---|
78 | real nox(nlayer) |
---|
79 | real no2x(nlayer) |
---|
80 | |
---|
81 | integer i,j,k,icol,indexint !indexes |
---|
82 | |
---|
83 | c variables for optical path calculation |
---|
84 | |
---|
85 | integer nz3 |
---|
86 | ! parameter (nz3=nz*2) |
---|
87 | |
---|
88 | integer jj |
---|
89 | real*8 esp(nlayer*2) |
---|
90 | real*8 ilayesp(nlayer*2) |
---|
91 | real*8 szalayesp(nlayer*2) |
---|
92 | integer nlayesp |
---|
93 | real*8 zmini |
---|
94 | real*8 depth |
---|
95 | real*8 espco2, espo2, espo3p, esph2, esph2o, esph2o2,espo3 |
---|
96 | real*8 espn2,espn,espno,espco,esph,espno2 |
---|
97 | real*8 rcmnz, rcmmini |
---|
98 | real*8 szadeg |
---|
99 | |
---|
100 | |
---|
101 | ! Tracer indexes in the thermospheric chemistry: |
---|
102 | !!! ATTENTION. These values have to be identical to those in chemthermos.F90 |
---|
103 | !!! If the values are changed there, the same has to be done here !!! |
---|
104 | integer,parameter :: i_co2 = 1 |
---|
105 | integer,parameter :: i_co = 2 |
---|
106 | integer,parameter :: i_o = 3 |
---|
107 | integer,parameter :: i_o1d = 4 |
---|
108 | integer,parameter :: i_o2 = 5 |
---|
109 | integer,parameter :: i_o3 = 6 |
---|
110 | integer,parameter :: i_h = 7 |
---|
111 | integer,parameter :: i_h2 = 8 |
---|
112 | integer,parameter :: i_oh = 9 |
---|
113 | integer,parameter :: i_ho2 = 10 |
---|
114 | integer,parameter :: i_h2o2 = 11 |
---|
115 | integer,parameter :: i_h2o = 12 |
---|
116 | integer,parameter :: i_n = 13 |
---|
117 | integer,parameter :: i_n2d = 14 |
---|
118 | integer,parameter :: i_no = 15 |
---|
119 | integer,parameter :: i_no2 = 16 |
---|
120 | integer,parameter :: i_n2 = 17 |
---|
121 | ! integer,parameter :: i_co2=1 |
---|
122 | ! integer,parameter :: i_o2=2 |
---|
123 | ! integer,parameter :: i_o=3 |
---|
124 | ! integer,parameter :: i_co=4 |
---|
125 | ! integer,parameter :: i_h=5 |
---|
126 | ! integer,parameter :: i_h2=8 |
---|
127 | ! integer,parameter :: i_h2o=9 |
---|
128 | ! integer,parameter :: i_h2o2=10 |
---|
129 | ! integer,parameter :: i_o3=12 |
---|
130 | ! integer,parameter :: i_n2=13 |
---|
131 | ! integer,parameter :: i_n=14 |
---|
132 | ! integer,parameter :: i_no=15 |
---|
133 | ! integer,parameter :: i_no2=17 |
---|
134 | |
---|
135 | |
---|
136 | c*************************PROGRAM STARTS******************************* |
---|
137 | |
---|
138 | nz3 = nlayer*2 |
---|
139 | do i=1,nlayer |
---|
140 | xx = ( radio + iz(i) ) * 1.e5 |
---|
141 | grav(i) = gg * masa /(xx**2) |
---|
142 | end do |
---|
143 | |
---|
144 | !Scale heights |
---|
145 | xx = kboltzman * tx(nlayer) * n_avog / grav(nlayer) |
---|
146 | Ho3p = xx / mmol(igcm_o) |
---|
147 | Hco2 = xx / mmol(igcm_co2) |
---|
148 | Ho2 = xx / mmol(igcm_o2) |
---|
149 | Hh2 = xx / mmol(igcm_h2) |
---|
150 | Hh2o = xx / mmol(igcm_h2o_vap) |
---|
151 | Hh2o2 = xx / mmol(igcm_h2o2) |
---|
152 | Hco = xx / mmol(igcm_co) |
---|
153 | Hh = xx / mmol(igcm_h) |
---|
154 | !Only if O3 chem. required |
---|
155 | if(chemthermod.ge.1) |
---|
156 | $ Ho3 = xx / mmol(igcm_o3) |
---|
157 | !Only if N or ion chem. |
---|
158 | if(chemthermod.ge.2) then |
---|
159 | Hn2 = xx / mmol(igcm_n2) |
---|
160 | Hn = xx / mmol(igcm_n) |
---|
161 | Hno = xx / mmol(igcm_no) |
---|
162 | Hno2 = xx / mmol(igcm_no2) |
---|
163 | endif |
---|
164 | ! first loop in altitude : initialisation |
---|
165 | do i=nlayer,1,-1 |
---|
166 | !Column initialisation |
---|
167 | co2colx(i) = 0. |
---|
168 | o2colx(i) = 0. |
---|
169 | o3pcolx(i) = 0. |
---|
170 | h2colx(i) = 0. |
---|
171 | h2ocolx(i) = 0. |
---|
172 | h2o2colx(i) = 0. |
---|
173 | o3colx(i) = 0. |
---|
174 | n2colx(i) = 0. |
---|
175 | ncolx(i) = 0. |
---|
176 | nocolx(i) = 0. |
---|
177 | cocolx(i) = 0. |
---|
178 | hcolx(i) = 0. |
---|
179 | no2colx(i) = 0. |
---|
180 | !Densities |
---|
181 | co2x(i) = rm(i,i_co2) |
---|
182 | o2x(i) = rm(i,i_o2) |
---|
183 | o3px(i) = rm(i,i_o) |
---|
184 | h2x(i) = rm(i,i_h2) |
---|
185 | h2ox(i) = rm(i,i_h2o) |
---|
186 | h2o2x(i) = rm(i,i_h2o2) |
---|
187 | cox(i) = rm(i,i_co) |
---|
188 | hx(i) = rm(i,i_h) |
---|
189 | !Only if O3 chem. required |
---|
190 | if(chemthermod.ge.1) |
---|
191 | $ o3x(i) = rm(i,i_o3) |
---|
192 | !Only if Nitrogen of ion chemistry requested |
---|
193 | if(chemthermod.ge.2) then |
---|
194 | n2x(i) = rm(i,i_n2) |
---|
195 | nx(i) = rm(i,i_n) |
---|
196 | nox(i) = rm(i,i_no) |
---|
197 | no2x(i) = rm(i,i_no2) |
---|
198 | endif |
---|
199 | enddo |
---|
200 | ! second loop in altitude : column calculations |
---|
201 | do i=nlayer,1,-1 |
---|
202 | !Routine to calculate the geometrical length of each layer |
---|
203 | call espesor_optico_A(ig,i,nlayer,zenit,iz(i),nz3,iz,esp, |
---|
204 | $ ilayesp,szalayesp,nlayesp, zmini) |
---|
205 | if(ilayesp(nlayesp).eq.-1) then |
---|
206 | co2colx(i)=1.e25 |
---|
207 | o2colx(i)=1.e25 |
---|
208 | o3pcolx(i)=1.e25 |
---|
209 | h2colx(i)=1.e25 |
---|
210 | h2ocolx(i)=1.e25 |
---|
211 | h2o2colx(i)=1.e25 |
---|
212 | o3colx(i)=1.e25 |
---|
213 | n2colx(i)=1.e25 |
---|
214 | ncolx(i)=1.e25 |
---|
215 | nocolx(i)=1.e25 |
---|
216 | cocolx(i)=1.e25 |
---|
217 | hcolx(i)=1.e25 |
---|
218 | no2colx(i)=1.e25 |
---|
219 | else |
---|
220 | rcmnz = ( radio + iz(nlayer) ) * 1.e5 |
---|
221 | rcmmini = ( radio + zmini ) * 1.e5 |
---|
222 | !Column calculation taking into account the geometrical depth |
---|
223 | !calculated before |
---|
224 | do j=1,nlayesp |
---|
225 | jj=ilayesp(j) |
---|
226 | !Top layer |
---|
227 | if(jj.eq.nlayer) then |
---|
228 | if(zenit.le.60.) then |
---|
229 | o3pcolx(i)=o3pcolx(i)+o3px(nlayer)*Ho3p*esp(j) |
---|
230 | $ *1.e-5 |
---|
231 | co2colx(i)=co2colx(i)+co2x(nlayer)*Hco2*esp(j) |
---|
232 | $ *1.e-5 |
---|
233 | h2o2colx(i)=h2o2colx(i)+ |
---|
234 | $ h2o2x(nlayer)*Hh2o2*esp(j)*1.e-5 |
---|
235 | o2colx(i)=o2colx(i)+o2x(nlayer)*Ho2*esp(j) |
---|
236 | $ *1.e-5 |
---|
237 | h2colx(i)=h2colx(i)+h2x(nlayer)*Hh2*esp(j) |
---|
238 | $ *1.e-5 |
---|
239 | h2ocolx(i)=h2ocolx(i)+h2ox(nlayer)*Hh2o*esp(j) |
---|
240 | $ *1.e-5 |
---|
241 | cocolx(i)=cocolx(i)+cox(nlayer)*Hco*esp(j) |
---|
242 | $ *1.e-5 |
---|
243 | hcolx(i)=hcolx(i)+hx(nlayer)*Hh*esp(j) |
---|
244 | $ *1.e-5 |
---|
245 | !Only if O3 chemistry required |
---|
246 | if(chemthermod.ge.1) o3colx(i)= |
---|
247 | $ o3colx(i)+o3x(nlayer)*Ho3*esp(j) |
---|
248 | $ *1.e-5 |
---|
249 | !Only if N or ion chemistry requested |
---|
250 | if(chemthermod.ge.2) then |
---|
251 | n2colx(i)=n2colx(i)+n2x(nlayer)*Hn2*esp(j) |
---|
252 | $ *1.e-5 |
---|
253 | ncolx(i)=ncolx(i)+nx(nlayer)*Hn*esp(j) |
---|
254 | $ *1.e-5 |
---|
255 | nocolx(i)=nocolx(i)+nox(nlayer)*Hno*esp(j) |
---|
256 | $ *1.e-5 |
---|
257 | no2colx(i)=no2colx(i)+no2x(nlayer)*Hno2*esp(j) |
---|
258 | $ *1.e-5 |
---|
259 | endif |
---|
260 | else if(zenit.gt.60.) then |
---|
261 | espco2 =sqrt((rcmnz+Hco2)**2 -rcmmini**2) - esp(j) |
---|
262 | espo2 = sqrt((rcmnz+Ho2)**2 -rcmmini**2) - esp(j) |
---|
263 | espo3p = sqrt((rcmnz+Ho3p)**2 -rcmmini**2)- esp(j) |
---|
264 | esph2 = sqrt((rcmnz+Hh2)**2 -rcmmini**2) - esp(j) |
---|
265 | esph2o = sqrt((rcmnz+Hh2o)**2 -rcmmini**2)- esp(j) |
---|
266 | esph2o2= sqrt((rcmnz+Hh2o2)**2-rcmmini**2)- esp(j) |
---|
267 | espco = sqrt((rcmnz+Hco)**2 -rcmmini**2) - esp(j) |
---|
268 | esph = sqrt((rcmnz+Hh)**2 -rcmmini**2) - esp(j) |
---|
269 | !Only if O3 chemistry required |
---|
270 | if(chemthermod.ge.1) |
---|
271 | $ espo3=sqrt((rcmnz+Ho3)**2-rcmmini**2)-esp(j) |
---|
272 | !Only if N or ion chemistry requested |
---|
273 | if(chemthermod.ge.2) then |
---|
274 | espn2 =sqrt((rcmnz+Hn2)**2-rcmmini**2)-esp(j) |
---|
275 | espn =sqrt((rcmnz+Hn)**2-rcmmini**2) - esp(j) |
---|
276 | espno =sqrt((rcmnz+Hno)**2-rcmmini**2) - esp(j) |
---|
277 | espno2=sqrt((rcmnz+Hno2)**2-rcmmini**2)- esp(j) |
---|
278 | endif |
---|
279 | co2colx(i) = co2colx(i) + espco2*co2x(nlayer) |
---|
280 | o2colx(i) = o2colx(i) + espo2*o2x(nlayer) |
---|
281 | o3pcolx(i) = o3pcolx(i) + espo3p*o3px(nlayer) |
---|
282 | h2colx(i) = h2colx(i) + esph2*h2x(nlayer) |
---|
283 | h2ocolx(i) = h2ocolx(i) + esph2o*h2ox(nlayer) |
---|
284 | h2o2colx(i)= h2o2colx(i)+ esph2o2*h2o2x(nlayer) |
---|
285 | cocolx(i) = cocolx(i) + espco*cox(nlayer) |
---|
286 | hcolx(i) = hcolx(i) + esph*hx(nlayer) |
---|
287 | !Only if O3 chemistry required |
---|
288 | if(chemthermod.ge.1) |
---|
289 | $ o3colx(i) = o3colx(i) + espo3*o3x(nlayer) |
---|
290 | !Only if N or ion chemistry requested |
---|
291 | if(chemthermod.ge.2) then |
---|
292 | n2colx(i) = n2colx(i) + espn2*n2x(nlayer) |
---|
293 | ncolx(i) = ncolx(i) + espn*nx(nlayer) |
---|
294 | nocolx(i) = nocolx(i) + espno*nox(nlayer) |
---|
295 | no2colx(i) = no2colx(i) + espno2*no2x(nlayer) |
---|
296 | endif |
---|
297 | endif !Of if zenit.lt.60 |
---|
298 | !Other layers |
---|
299 | else |
---|
300 | co2colx(i) = co2colx(i) + |
---|
301 | $ esp(j) * (co2x(jj)+co2x(jj+1)) / 2. |
---|
302 | o2colx(i) = o2colx(i) + |
---|
303 | $ esp(j) * (o2x(jj)+o2x(jj+1)) / 2. |
---|
304 | o3pcolx(i) = o3pcolx(i) + |
---|
305 | $ esp(j) * (o3px(jj)+o3px(jj+1)) / 2. |
---|
306 | h2colx(i) = h2colx(i) + |
---|
307 | $ esp(j) * (h2x(jj)+h2x(jj+1)) / 2. |
---|
308 | h2ocolx(i) = h2ocolx(i) + |
---|
309 | $ esp(j) * (h2ox(jj)+h2ox(jj+1)) / 2. |
---|
310 | h2o2colx(i) = h2o2colx(i) + |
---|
311 | $ esp(j) * (h2o2x(jj)+h2o2x(jj+1)) / 2. |
---|
312 | cocolx(i) = cocolx(i) + |
---|
313 | $ esp(j) * (cox(jj)+cox(jj+1)) / 2. |
---|
314 | hcolx(i) = hcolx(i) + |
---|
315 | $ esp(j) * (hx(jj)+hx(jj+1)) / 2. |
---|
316 | !Only if O3 chemistry required |
---|
317 | if(chemthermod.ge.1) |
---|
318 | $ o3colx(i) = o3colx(i) + |
---|
319 | $ esp(j) * (o3x(jj)+o3x(jj+1)) / 2. |
---|
320 | !Only if N or ion chemistry requested |
---|
321 | if(chemthermod.ge.2) then |
---|
322 | n2colx(i) = n2colx(i) + |
---|
323 | $ esp(j) * (n2x(jj)+n2x(jj+1)) / 2. |
---|
324 | ncolx(i) = ncolx(i) + |
---|
325 | $ esp(j) * (nx(jj)+nx(jj+1)) / 2. |
---|
326 | nocolx(i) = nocolx(i) + |
---|
327 | $ esp(j) * (nox(jj)+nox(jj+1)) / 2. |
---|
328 | no2colx(i) = no2colx(i) + |
---|
329 | $ esp(j) * (no2x(jj)+no2x(jj+1)) / 2. |
---|
330 | endif |
---|
331 | endif !Of if jj.eq.nlayer |
---|
332 | end do !Of do j=1,nlayesp |
---|
333 | endif !Of ilayesp(nlayesp).eq.-1 |
---|
334 | enddo !Of do i=nlayer,1,-1 |
---|
335 | |
---|
336 | |
---|
337 | end subroutine column |
---|
338 | |
---|
339 | |
---|
340 | c********************************************************************** |
---|
341 | c********************************************************************** |
---|
342 | |
---|
343 | subroutine interfast(wm,wp,nm,p,nlayer,pin,nl,limdown,limup) |
---|
344 | C |
---|
345 | C subroutine to perform linear interpolation in pressure from 1D profile |
---|
346 | C escin(nl) sampled on pressure grid pin(nl) to profile |
---|
347 | C escout(nlayer) on pressure grid p(nlayer). |
---|
348 | C |
---|
349 | real*8,intent(out) :: wm(nlayer),wp(nlayer) ! interpolation weights |
---|
350 | integer,intent(out) :: nm(nlayer) ! index of nearest point |
---|
351 | real*8,intent(in) :: pin(nl),p(nlayer) |
---|
352 | real*8,intent(in) :: limup,limdown |
---|
353 | integer,intent(in) :: nl,nlayer |
---|
354 | integer :: n1,n,np,nini |
---|
355 | logical,parameter :: extra_sanity_checks=.false. |
---|
356 | |
---|
357 | ! Added sanity check: is input p(:) indeed monotonically increasing? |
---|
358 | if (extra_sanity_checks) then |
---|
359 | do nini=1,nlayer-1 |
---|
360 | if (p(nini).gt.p(nini+1)) then |
---|
361 | write(*,*) "possible interfast issue, nini=",nini |
---|
362 | write(*,*) "p(nini)=",p(nini),"> p(nini+1)=",p(nini+1) |
---|
363 | endif |
---|
364 | enddo |
---|
365 | endif ! of if (extra_sanity_checks) |
---|
366 | |
---|
367 | nini=1 |
---|
368 | do n1=1,nlayer |
---|
369 | if(p(n1) .gt. limup .or. p(n1) .lt. limdown) then |
---|
370 | wm(n1) = 0.d0 |
---|
371 | wp(n1) = 0.d0 |
---|
372 | else |
---|
373 | do n = nini,nl-1 |
---|
374 | if (p(n1).ge.pin(n).and.p(n1).le.pin(n+1)) then |
---|
375 | nm(n1)=n |
---|
376 | np=n+1 |
---|
377 | wm(n1)=abs(pin(n)-p(n1))/(pin(np)-pin(n)) |
---|
378 | wp(n1)=1.d0 - wm(n1) |
---|
379 | nini = n |
---|
380 | exit |
---|
381 | endif |
---|
382 | enddo |
---|
383 | endif |
---|
384 | enddo |
---|
385 | |
---|
386 | ! Added sanity check: does nm(:) indeed contain values |
---|
387 | ! between 1 and nl-1 ? |
---|
388 | if (extra_sanity_checks) then |
---|
389 | if ((minval(nm)<1).or.(maxval(nm)>nl-1)) then |
---|
390 | write(*,*) "interfast issue nm(:) contains incoherent values" |
---|
391 | write(*,*) " nm(:) values should be between 1 and ",nl-1 |
---|
392 | write(*,*) " but nm(:)=",nm(:) |
---|
393 | endif |
---|
394 | endif ! of if (extra_sanity_checks) |
---|
395 | |
---|
396 | end subroutine interfast |
---|
397 | |
---|
398 | |
---|
399 | c********************************************************************** |
---|
400 | c********************************************************************** |
---|
401 | |
---|
402 | subroutine espesor_optico_A (ig,capa,nlayer, szadeg,z, |
---|
403 | @ nz3,iz,esp,ilayesp,szalayesp,nlayesp, zmini) |
---|
404 | |
---|
405 | c fgg nov 03 Adaptation to Martian model |
---|
406 | c malv jul 03 Corrected z grid. Split in alt & frec codes |
---|
407 | c fgg feb 03 first version |
---|
408 | ************************************************************************* |
---|
409 | |
---|
410 | use param_v4_h, only: radio |
---|
411 | implicit none |
---|
412 | |
---|
413 | c arguments |
---|
414 | |
---|
415 | real szadeg ! I. SZA [rad] |
---|
416 | real z ! I. altitude of interest [km] |
---|
417 | integer nz3,ig,nlayer ! I. dimension of esp, ylayesp, etc... |
---|
418 | ! (=2*nlayer= max# of layers in ray path) |
---|
419 | real iz(nlayer+1) ! I. Altitude of each layer |
---|
420 | real*8 esp(nz3) ! O. layer widths after geometrically |
---|
421 | ! amplified; in [cm] except at TOA |
---|
422 | ! where an auxiliary value is used |
---|
423 | real*8 ilayesp(nz3) ! O. Indexes of layers along ray path |
---|
424 | real*8 szalayesp(nz3) ! O. SZA [deg] " " " |
---|
425 | integer nlayesp |
---|
426 | ! real*8 nlayesp ! O. # layers along ray path at this z |
---|
427 | real*8 zmini ! O. Minimum altitud of ray path [km] |
---|
428 | |
---|
429 | |
---|
430 | c local variables and constants |
---|
431 | |
---|
432 | integer j,i,capa |
---|
433 | integer jmin ! index of min.altitude along ray path |
---|
434 | real*8 szarad ! SZA [deg] |
---|
435 | real*8 zz |
---|
436 | real*8 diz(nlayer+1) |
---|
437 | real*8 rkmnz ! distance TOA to center of Planet [km] |
---|
438 | real*8 rkmmini ! distance zmini to center of P [km] |
---|
439 | real*8 rkmj ! intermediate distance to C of P [km] |
---|
440 | c external function |
---|
441 | c external grid_R8 ! Returns index of layer containing the altitude |
---|
442 | c ! of interest, z; for example, if |
---|
443 | c ! zkm(i)=z or zkm(i)<z<zkm(i+1) => grid(z)=i |
---|
444 | c integer grid_R8 |
---|
445 | |
---|
446 | ************************************************************************* |
---|
447 | szarad = dble(szadeg)*3.141592d0/180.d0 |
---|
448 | zz=dble(z) |
---|
449 | do i=1,nlayer |
---|
450 | diz(i)=dble(iz(i)) |
---|
451 | enddo |
---|
452 | do j=1,nz3 |
---|
453 | esp(j) = 0.d0 |
---|
454 | szalayesp(j) = 777.d0 |
---|
455 | ilayesp(j) = 0 |
---|
456 | enddo |
---|
457 | nlayesp = 0 |
---|
458 | |
---|
459 | ! First case: szadeg<60 |
---|
460 | ! The optical thickness will be given by 1/cos(sza) |
---|
461 | ! We deal with 2 different regions: |
---|
462 | ! 1: First, all layers between z and ztop ("upper part of ray") |
---|
463 | ! 2: Second, the layer at ztop |
---|
464 | if(szadeg.lt.60.d0) then |
---|
465 | |
---|
466 | zmini = zz |
---|
467 | if(abs(zz-diz(nlayer)).lt.1.d-3) goto 1357 |
---|
468 | ! 1st Zone: Upper part of ray |
---|
469 | ! |
---|
470 | do j=grid_R8(zz,diz,nlayer),nlayer-1 |
---|
471 | nlayesp = nlayesp + 1 |
---|
472 | ilayesp(nlayesp) = j |
---|
473 | esp(nlayesp) = (diz(j+1)-diz(j)) / cos(szarad) ! [km] |
---|
474 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
475 | szalayesp(nlayesp) = szadeg |
---|
476 | end do |
---|
477 | |
---|
478 | ! |
---|
479 | ! 2nd Zone: Top layer |
---|
480 | 1357 continue |
---|
481 | nlayesp = nlayesp + 1 |
---|
482 | ilayesp(nlayesp) = nlayer |
---|
483 | esp(nlayesp) = 1.d0 / cos(szarad) ! aux. non-dimens. factor |
---|
484 | szalayesp(nlayesp) = szadeg |
---|
485 | |
---|
486 | |
---|
487 | ! Second case: 60 < szadeg < 90 |
---|
488 | ! The optical thickness is evaluated. |
---|
489 | ! (the magnitude of the effect of not using cos(sza) is 3.e-5 |
---|
490 | ! for z=60km & sza=30, and 5e-4 for z=60km & sza=60, approximately) |
---|
491 | ! We deal with 2 different regions: |
---|
492 | ! 1: First, all layers between z and ztop ("upper part of ray") |
---|
493 | ! 2: Second, the layer at ztop ("uppermost layer") |
---|
494 | else if(szadeg.le.90.d0.and.szadeg.ge.60.d0) then |
---|
495 | |
---|
496 | zmini=(radio+zz)*sin(szarad)-radio |
---|
497 | rkmmini = radio + zmini |
---|
498 | |
---|
499 | if(abs(zz-diz(nlayer)).lt.1.d-4) goto 1470 |
---|
500 | |
---|
501 | ! 1st Zone: Upper part of ray |
---|
502 | ! |
---|
503 | do j=grid_R8(zz,diz,nlayer),nlayer-1 |
---|
504 | nlayesp = nlayesp + 1 |
---|
505 | ilayesp(nlayesp) = j |
---|
506 | esp(nlayesp) = |
---|
507 | & sqrt( (radio+diz(j+1))**2 - rkmmini**2 ) - |
---|
508 | & sqrt( (radio+diz(j))**2 - rkmmini**2 ) ! [km] |
---|
509 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
510 | rkmj = radio+diz(j) |
---|
511 | szalayesp(nlayesp) = asin( rkmmini/rkmj ) ! [rad] |
---|
512 | szalayesp(nlayesp) = szalayesp(nlayesp) * 180.d0/3.141592 ! [deg] |
---|
513 | end do |
---|
514 | 1470 continue |
---|
515 | ! 2nd Zone: Uppermost layer of ray. |
---|
516 | ! |
---|
517 | nlayesp = nlayesp + 1 |
---|
518 | ilayesp(nlayesp) = nlayer |
---|
519 | rkmnz = radio+diz(nlayer) |
---|
520 | esp(nlayesp) = sqrt( rkmnz**2 - rkmmini**2 ) ! aux.factor[km] |
---|
521 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! aux.f. [cm] |
---|
522 | szalayesp(nlayesp) = asin( rkmmini/rkmnz ) ! [rad] |
---|
523 | szalayesp(nlayesp) = szalayesp(nlayesp) * 180.d0/3.141592! [deg] |
---|
524 | |
---|
525 | |
---|
526 | ! Third case: szadeg > 90 |
---|
527 | ! The optical thickness is evaluated. |
---|
528 | ! We deal with 5 different regions: |
---|
529 | ! 1: all layers between z and ztop ("upper part of ray") |
---|
530 | ! 2: the layer at ztop ("uppermost layer") |
---|
531 | ! 3: the lowest layer, at zmini |
---|
532 | ! 4: the layers increasing from zmini to z (here SZA<90) |
---|
533 | ! 5: the layers decreasing from z to zmini (here SZA>90) |
---|
534 | else if(szadeg.gt.90.d0) then |
---|
535 | |
---|
536 | zmini=(radio+zz)*sin(szarad)-radio |
---|
537 | !zmini should be lower than zz, as SZA<90. However, in situations |
---|
538 | !where SZA is very close to 90, rounding errors can make zmini |
---|
539 | !slightly higher than zz, causing problems in the determination |
---|
540 | !of the jmin index. A correction is implemented in the determination |
---|
541 | !of jmin, some lines below |
---|
542 | rkmmini = radio + zmini |
---|
543 | |
---|
544 | if(zmini.lt.diz(1)) then ! Can see the sun? No => esp(j)=inft |
---|
545 | nlayesp = nlayesp + 1 |
---|
546 | ilayesp(nlayesp) = - 1 ! Value to mark "no sun on view" |
---|
547 | ! esp(nlayesp) = 1.e30 |
---|
548 | |
---|
549 | else |
---|
550 | jmin=grid_R8(zmini,diz,nlayer)+1 |
---|
551 | !Correction for possible rounding errors when SZA very close |
---|
552 | !to 90 degrees |
---|
553 | if(jmin.gt.grid_R8(zz,diz,nlayer)) then |
---|
554 | write(*,*)'jthermcalc warning: possible rounding error' |
---|
555 | write(*,*)'point,sza,layer:',ig,szadeg,capa |
---|
556 | jmin=grid_R8(zz,diz,nlayer) |
---|
557 | endif |
---|
558 | |
---|
559 | if(abs(zz-diz(nlayer)).lt.1.d-4) goto 9876 |
---|
560 | |
---|
561 | ! 1st Zone: Upper part of ray |
---|
562 | ! |
---|
563 | do j=grid_R8(zz,diz,nlayer),nlayer-1 |
---|
564 | nlayesp = nlayesp + 1 |
---|
565 | ilayesp(nlayesp) = j |
---|
566 | esp(nlayesp) = |
---|
567 | $ sqrt( (radio+diz(j+1))**2 - rkmmini**2 ) - |
---|
568 | $ sqrt( (radio+diz(j))**2 - rkmmini**2 ) ! [km] |
---|
569 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
570 | rkmj = radio+diz(j) |
---|
571 | szalayesp(nlayesp) = asin( rkmmini/rkmj ) ! [rad] |
---|
572 | szalayesp(nlayesp) = szalayesp(nlayesp) *180.d0/3.141592 ! [deg] |
---|
573 | end do |
---|
574 | |
---|
575 | 9876 continue |
---|
576 | ! 2nd Zone: Uppermost layer of ray. |
---|
577 | ! |
---|
578 | nlayesp = nlayesp + 1 |
---|
579 | ilayesp(nlayesp) = nlayer |
---|
580 | rkmnz = radio+diz(nlayer) |
---|
581 | esp(nlayesp) = sqrt( rkmnz**2 - rkmmini**2 ) !aux.factor[km] |
---|
582 | esp(nlayesp) = esp(nlayesp) * 1.d5 !aux.f.[cm] |
---|
583 | szalayesp(nlayesp) = asin( rkmmini/rkmnz ) ! [rad] |
---|
584 | szalayesp(nlayesp) = szalayesp(nlayesp) *180.d0/3.141592 ! [deg] |
---|
585 | |
---|
586 | ! 3er Zone: Lowestmost layer of ray |
---|
587 | ! |
---|
588 | if ( jmin .ge. 2 ) then ! If above the planet's surface |
---|
589 | j=jmin-1 |
---|
590 | nlayesp = nlayesp + 1 |
---|
591 | ilayesp(nlayesp) = j |
---|
592 | esp(nlayesp) = 2. * |
---|
593 | $ sqrt( (radio+diz(j+1))**2 -rkmmini**2 ) ! [km] |
---|
594 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
595 | rkmj = radio+diz(j+1) |
---|
596 | szalayesp(nlayesp) = asin( rkmmini/rkmj ) ! [rad] |
---|
597 | szalayesp(nlayesp) = szalayesp(nlayesp) *180.d0/3.141592 ! [deg] |
---|
598 | endif |
---|
599 | |
---|
600 | ! 4th zone: Lower part of ray, increasing from zmin to z |
---|
601 | ! ( layers with SZA < 90 deg ) |
---|
602 | do j=jmin,grid_R8(zz,diz,nlayer)-1 |
---|
603 | nlayesp = nlayesp + 1 |
---|
604 | ilayesp(nlayesp) = j |
---|
605 | esp(nlayesp) = |
---|
606 | $ sqrt( (radio+diz(j+1))**2 - rkmmini**2 ) |
---|
607 | $ - sqrt( (radio+diz(j))**2 - rkmmini**2 ) ! [km] |
---|
608 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
609 | rkmj = radio+diz(j) |
---|
610 | szalayesp(nlayesp) = asin( rkmmini/rkmj ) ! [rad] |
---|
611 | szalayesp(nlayesp) = szalayesp(nlayesp) *180.d0/3.141592 ! [deg] |
---|
612 | end do |
---|
613 | |
---|
614 | ! 5th zone: Lower part of ray, decreasing from z to zmin |
---|
615 | ! ( layers with SZA > 90 deg ) |
---|
616 | do j=grid_R8(zz,diz,nlayer)-1, jmin, -1 |
---|
617 | nlayesp = nlayesp + 1 |
---|
618 | ilayesp(nlayesp) = j |
---|
619 | esp(nlayesp) = |
---|
620 | $ sqrt( (radio+diz(j+1))**2 - rkmmini**2 ) |
---|
621 | $ - sqrt( (radio+diz(j))**2 - rkmmini**2 ) ! [km] |
---|
622 | esp(nlayesp) = esp(nlayesp) * 1.d5 ! [cm] |
---|
623 | rkmj = radio+diz(j) |
---|
624 | szalayesp(nlayesp) = 3.141592 - asin( rkmmini/rkmj ) ! [rad] |
---|
625 | szalayesp(nlayesp) = szalayesp(nlayesp)*180.d0/3.141592 ! [deg] |
---|
626 | end do |
---|
627 | |
---|
628 | end if |
---|
629 | |
---|
630 | end if |
---|
631 | |
---|
632 | |
---|
633 | end subroutine espesor_optico_A |
---|
634 | |
---|
635 | |
---|
636 | |
---|
637 | c********************************************************************** |
---|
638 | c*********************************************************************** |
---|
639 | |
---|
640 | function grid_R8 (z, zgrid, nz) |
---|
641 | |
---|
642 | c Returns the index where z is located within vector zgrid |
---|
643 | c The vector zgrid must be monotonously increasing, otherwise program stops. |
---|
644 | c If z is outside zgrid limits, or zgrid dimension is nz<2, the program stops. |
---|
645 | c |
---|
646 | c FGG Aug-2004 Correct z.lt.zgrid(i) to .le. |
---|
647 | c MALV Jul-2003 |
---|
648 | c*********************************************************************** |
---|
649 | |
---|
650 | implicit none |
---|
651 | |
---|
652 | c Arguments |
---|
653 | integer nz |
---|
654 | real*8 z |
---|
655 | real*8 zgrid(nz) |
---|
656 | integer grid_R8 |
---|
657 | |
---|
658 | c Local |
---|
659 | integer i, nz1, nznew |
---|
660 | |
---|
661 | c*** CODE START |
---|
662 | |
---|
663 | if ( z .lt. zgrid(1) ) then |
---|
664 | write (*,*) ' GRID/ z outside bounds of zgrid ' |
---|
665 | write (*,*) ' z,zgrid(1),zgrid(nz) =', z,zgrid(1),zgrid(nz) |
---|
666 | z = zgrid(1) |
---|
667 | write(*,*) 'WARNING: error in grid_r8 (jthermcalc.F)' |
---|
668 | write(*,*) 'Please check values of z and zgrid above' |
---|
669 | endif |
---|
670 | if (z .gt. zgrid(nz) ) then |
---|
671 | write (*,*) ' GRID/ z outside bounds of zgrid ' |
---|
672 | write (*,*) ' z,zgrid(1),zgrid(nz) =', z,zgrid(1),zgrid(nz) |
---|
673 | z = zgrid(nz) |
---|
674 | write(*,*) 'WARNING: error in grid_r8 (jthermcalc.F)' |
---|
675 | write(*,*) 'Please check values of z and zgrid above' |
---|
676 | endif |
---|
677 | if ( nz .lt. 2 ) then |
---|
678 | write (*,*) ' GRID/ zgrid needs 2 points at least ! ' |
---|
679 | stop ' Serious error in GRID.F ' |
---|
680 | endif |
---|
681 | if ( zgrid(1) .ge. zgrid(nz) ) then |
---|
682 | write (*,*) ' GRID/ zgrid must increase with index' |
---|
683 | stop ' Serious error in GRID.F ' |
---|
684 | endif |
---|
685 | |
---|
686 | nz1 = 1 |
---|
687 | nznew = nz/2 |
---|
688 | if ( z .gt. zgrid(nznew) ) then |
---|
689 | nz1 = nznew |
---|
690 | nznew = nz |
---|
691 | endif |
---|
692 | do i=nz1+1,nznew |
---|
693 | if ( z. eq. zgrid(i) ) then |
---|
694 | grid_R8=i |
---|
695 | return |
---|
696 | elseif ( z .le. zgrid(i) ) then |
---|
697 | grid_R8 = i-1 |
---|
698 | return |
---|
699 | endif |
---|
700 | enddo |
---|
701 | grid_R8 = nz |
---|
702 | |
---|
703 | |
---|
704 | |
---|
705 | end function grid_R8 |
---|
706 | |
---|
707 | |
---|
708 | |
---|
709 | !c*************************************************** |
---|
710 | !c*************************************************** |
---|
711 | |
---|
712 | subroutine flujo(date) |
---|
713 | |
---|
714 | |
---|
715 | !c fgg nov 2002 first version |
---|
716 | !c*************************************************** |
---|
717 | |
---|
718 | use comsaison_h, only: dist_sol |
---|
719 | use param_v4_h, only: ninter, |
---|
720 | . fluxtop, ct1, ct2, p1, p2 |
---|
721 | implicit none |
---|
722 | |
---|
723 | |
---|
724 | ! common variables and constants |
---|
725 | include "callkeys.h" |
---|
726 | |
---|
727 | |
---|
728 | ! Arguments |
---|
729 | |
---|
730 | real date |
---|
731 | |
---|
732 | |
---|
733 | ! Local variable and constants |
---|
734 | |
---|
735 | integer i |
---|
736 | integer inter |
---|
737 | real nada |
---|
738 | |
---|
739 | !c************************************************* |
---|
740 | |
---|
741 | if(date.lt.1985.) date=1985. |
---|
742 | if(date.gt.2001.) date=2001. |
---|
743 | |
---|
744 | do i=1,ninter |
---|
745 | fluxtop(i)=1. |
---|
746 | !Variation of solar flux with 11 years solar cycle |
---|
747 | !For more details, see Gonzalez-Galindo et al. 2005 |
---|
748 | !To be improved in next versions |
---|
749 | if(i.le.24.and.solvarmod.eq.0) then |
---|
750 | fluxtop(i)=(((ct1(i)+p1(i)*date)/2.) |
---|
751 | $ *sin(2.*3.1416/11.*(date-1985.-3.1416)) |
---|
752 | $ +(ct2(i)+p2(i)*date)+1.)*fluxtop(i) |
---|
753 | end if |
---|
754 | fluxtop(i)=fluxtop(i)*(1.52/dist_sol)**2 |
---|
755 | end do |
---|
756 | |
---|
757 | end subroutine flujo |
---|
758 | |
---|
759 | end module jthermcalc_util |
---|