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