1 | subroutine interpolateH2He(wn,temp,presH2,presHe,abcoef,firstcall,ind) |
---|
2 | |
---|
3 | !================================================================== |
---|
4 | ! |
---|
5 | ! Purpose |
---|
6 | ! ------- |
---|
7 | ! Calculates the H2-He CIA opacity, using a lookup table from |
---|
8 | ! HITRAN (2011) |
---|
9 | ! |
---|
10 | ! Authors |
---|
11 | ! ------- |
---|
12 | ! R. Wordsworth (2011) |
---|
13 | ! |
---|
14 | !================================================================== |
---|
15 | |
---|
16 | use datafile_mod, only: datadir |
---|
17 | |
---|
18 | implicit none |
---|
19 | |
---|
20 | ! input |
---|
21 | double precision wn ! wavenumber (cm^-1) |
---|
22 | double precision temp ! temperature (Kelvin) |
---|
23 | double precision presH2 ! H2 partial pressure (Pascals) |
---|
24 | double precision presHe ! He partial pressure (Pascals) |
---|
25 | |
---|
26 | ! output |
---|
27 | double precision abcoef ! absorption coefficient (m^-1) |
---|
28 | |
---|
29 | integer nS,nT |
---|
30 | parameter(nS=2428) |
---|
31 | parameter(nT=10) |
---|
32 | |
---|
33 | double precision, parameter :: losch = 2.6867774e19 |
---|
34 | ! Loschmit's number (molecule cm^-3 at STP) |
---|
35 | ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2 |
---|
36 | ! see Richard et al. 2011, JQSRT for details |
---|
37 | |
---|
38 | double precision amagatH2 |
---|
39 | double precision amagatHe |
---|
40 | double precision wn_arr(nS) |
---|
41 | double precision temp_arr(nT) |
---|
42 | double precision abs_arr(nS,nT) |
---|
43 | |
---|
44 | integer k,iT |
---|
45 | logical firstcall |
---|
46 | |
---|
47 | save wn_arr, temp_arr, abs_arr |
---|
48 | |
---|
49 | character*100 dt_file |
---|
50 | integer strlen,ios |
---|
51 | |
---|
52 | character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)" |
---|
53 | |
---|
54 | character*20 bleh |
---|
55 | double precision blah, Ttemp |
---|
56 | integer nres |
---|
57 | |
---|
58 | integer ind |
---|
59 | |
---|
60 | if(temp.gt.400)then |
---|
61 | print*,'Your temperatures are too high for this H2-He CIA dataset.' |
---|
62 | print*,'Please run mixed H2-He atmospheres below T = 400 K.' |
---|
63 | stop |
---|
64 | endif |
---|
65 | |
---|
66 | amagatH2 = (273.15/temp)*(presH2/101325.0) |
---|
67 | amagatHe = (273.15/temp)*(presHe/101325.0) |
---|
68 | |
---|
69 | if(firstcall)then ! called by sugas_corrk only |
---|
70 | print*,'----------------------------------------------------' |
---|
71 | print*,'Initialising H2-He continuum from HITRAN database...' |
---|
72 | |
---|
73 | ! 1.1 Open the ASCII files |
---|
74 | dt_file=TRIM(datadir)//'/continuum_data/H2-He_norm_2011.cia' |
---|
75 | |
---|
76 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
77 | if (ios.ne.0) then ! file not found |
---|
78 | write(*,*) 'Error from interpolateH2He' |
---|
79 | write(*,*) 'Data file ',trim(dt_file),' not found.' |
---|
80 | write(*,*) 'Check that your path to datagcm:',trim(datadir) |
---|
81 | write(*,*) 'is correct. You can change it in callphys.def with:' |
---|
82 | write(*,*) 'datadir = /absolute/path/to/datagcm' |
---|
83 | write(*,*) 'Also check that the continuum data continuum_data/H2-He_norm_2011.cia is there.' |
---|
84 | call abort |
---|
85 | else |
---|
86 | |
---|
87 | do iT=1,nT |
---|
88 | |
---|
89 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
90 | if(nS.ne.nres)then |
---|
91 | print*,'Resolution given in file: ',trim(dt_file) |
---|
92 | print*,'is ',nres,', which does not match nS.' |
---|
93 | print*,'Please adjust nS value in interpolateH2He.F90' |
---|
94 | stop |
---|
95 | endif |
---|
96 | temp_arr(iT)=Ttemp |
---|
97 | |
---|
98 | do k=1,nS |
---|
99 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
100 | end do |
---|
101 | |
---|
102 | end do |
---|
103 | |
---|
104 | endif |
---|
105 | close(33) |
---|
106 | |
---|
107 | print*,'interpolateH2He: At wavenumber ',wn,' cm^-1' |
---|
108 | print*,' temperature ',temp,' K' |
---|
109 | print*,' H2 partial pressure ',presH2,' Pa' |
---|
110 | print*,' and He partial pressure ',presHe,' Pa' |
---|
111 | |
---|
112 | endif |
---|
113 | |
---|
114 | call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind) |
---|
115 | |
---|
116 | !print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
117 | !print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
118 | |
---|
119 | abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1 |
---|
120 | |
---|
121 | !print*,'We have ',amagatH2,' amagats of H2' |
---|
122 | !print*,'and ',amagatHe,' amagats of He' |
---|
123 | !print*,'So the absorption is ',abcoef,' m^-1' |
---|
124 | |
---|
125 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
126 | ! however our bands are normally thin, so this is no big deal. |
---|
127 | |
---|
128 | return |
---|
129 | end subroutine interpolateH2He |
---|