[253] | 1 | subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall) |
---|
| 2 | |
---|
| 3 | !================================================================== |
---|
| 4 | ! |
---|
| 5 | ! Purpose |
---|
| 6 | ! ------- |
---|
| 7 | ! Calculates the H2-H2 CIA opacity, using a lookup table from |
---|
[716] | 8 | ! HITRAN (2011) |
---|
[253] | 9 | ! |
---|
| 10 | ! Authors |
---|
| 11 | ! ------- |
---|
[716] | 12 | ! R. Wordsworth (2011) |
---|
[253] | 13 | ! |
---|
| 14 | !================================================================== |
---|
| 15 | |
---|
[374] | 16 | use datafile_mod, only: datadir |
---|
| 17 | implicit none |
---|
[253] | 18 | |
---|
| 19 | ! input |
---|
| 20 | double precision wn ! wavenumber (cm^-1) |
---|
| 21 | double precision temp ! temperature (Kelvin) |
---|
| 22 | double precision pres ! pressure (Pascals) |
---|
| 23 | |
---|
| 24 | ! output |
---|
| 25 | double precision abcoef ! absorption coefficient (m^-1) |
---|
| 26 | |
---|
| 27 | integer nS,nT |
---|
[716] | 28 | parameter(nS=2428) |
---|
| 29 | parameter(nT=10) |
---|
[253] | 30 | |
---|
[716] | 31 | double precision, parameter :: losch = 2.6867774e19 |
---|
| 32 | ! Loschmit's number (molecule cm^-3 at STP) |
---|
| 33 | ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2 |
---|
| 34 | ! see Richard et al. 2011, JQSRT for details |
---|
| 35 | |
---|
[253] | 36 | double precision amagat |
---|
| 37 | double precision wn_arr(nS) |
---|
| 38 | double precision temp_arr(nT) |
---|
| 39 | double precision abs_arr(nS,nT) |
---|
| 40 | |
---|
[716] | 41 | integer k,iT |
---|
[253] | 42 | logical firstcall |
---|
| 43 | |
---|
| 44 | save wn_arr, temp_arr, abs_arr |
---|
| 45 | |
---|
| 46 | character*100 dt_file |
---|
| 47 | integer strlen,ios |
---|
| 48 | |
---|
[716] | 49 | character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)" |
---|
[253] | 50 | |
---|
[716] | 51 | character*20 bleh |
---|
| 52 | double precision blah, Ttemp |
---|
| 53 | integer nres |
---|
[253] | 54 | |
---|
[716] | 55 | if(temp.gt.400)then |
---|
| 56 | print*,'Your temperatures are too high for this H2-H2 CIA dataset. If you ' |
---|
| 57 | print*,'really want to run simulations with hydrogen at T > 400 K, contact' |
---|
| 58 | print*,'Robin Wordsworth [rwordsworth@uchicago.edu].' |
---|
| 59 | stop |
---|
| 60 | endif |
---|
| 61 | |
---|
| 62 | amagat = (273.15/temp)*(pres/101325.0) |
---|
| 63 | |
---|
| 64 | if(firstcall)then ! called by sugas_corrk only |
---|
| 65 | print*,'----------------------------------------------------' |
---|
| 66 | print*,'Initialising H2-H2 continuum from HITRAN database...' |
---|
| 67 | |
---|
[253] | 68 | ! 1.1 Open the ASCII files |
---|
[716] | 69 | dt_file=TRIM(datadir)//'/continuum_data/H2-H2_norm_2011.cia' |
---|
[253] | 70 | |
---|
| 71 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
| 72 | if (ios.ne.0) then ! file not found |
---|
[374] | 73 | write(*,*) 'Error from interpolateH2H2' |
---|
| 74 | write(*,*) 'Data file ',trim(dt_file),' not found.' |
---|
[716] | 75 | write(*,*) 'Check that your path to datagcm:',trim(datadir) |
---|
| 76 | write(*,*) 'is correct. You can change it in callphys.def with:' |
---|
| 77 | write(*,*) 'datadir = /absolute/path/to/datagcm' |
---|
| 78 | write(*,*) 'Also check that the continuum data continuum_data/H2-H2_norm_2011.cia is there.' |
---|
[374] | 79 | call abort |
---|
[253] | 80 | else |
---|
[716] | 81 | |
---|
| 82 | do iT=1,nT |
---|
| 83 | |
---|
| 84 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
| 85 | if(nS.ne.nres)then |
---|
| 86 | print*,'Resolution given in file: ',trim(dt_file) |
---|
| 87 | print*,'is ',nres,', which does not match nS.' |
---|
| 88 | print*,'Please adjust nS value in interpolateH2H2.F90' |
---|
| 89 | stop |
---|
| 90 | endif |
---|
| 91 | temp_arr(iT)=Ttemp |
---|
| 92 | |
---|
| 93 | do k=1,nS |
---|
| 94 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
| 95 | end do |
---|
| 96 | |
---|
[253] | 97 | end do |
---|
[716] | 98 | |
---|
[253] | 99 | endif |
---|
| 100 | close(33) |
---|
| 101 | |
---|
[374] | 102 | print*,'interpolateH2H2: At wavenumber ',wn,' cm^-1' |
---|
[253] | 103 | print*,' temperature ',temp,' K' |
---|
| 104 | print*,' pressure ',pres,' Pa' |
---|
| 105 | |
---|
[305] | 106 | call bilinearH2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef) |
---|
[253] | 107 | |
---|
[716] | 108 | print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
| 109 | print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
[253] | 110 | |
---|
[716] | 111 | abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1 |
---|
[253] | 112 | |
---|
[716] | 113 | print*,'We have ',amagat,' amagats of H2' |
---|
[253] | 114 | print*,'So the absorption is ',abcoef,' m^-1' |
---|
| 115 | |
---|
[716] | 116 | !open(117,file='T_array.dat') |
---|
| 117 | !do iT=1,nT |
---|
| 118 | ! write(117,*), temp_arr(iT) |
---|
| 119 | !end do |
---|
| 120 | !close(117) |
---|
| 121 | |
---|
| 122 | !open(115,file='nu_array.dat') |
---|
| 123 | !do k=1,nS |
---|
| 124 | ! write(115,*), wn_arr(k) |
---|
| 125 | !end do |
---|
| 126 | !close(115) |
---|
| 127 | |
---|
| 128 | !open(113,file='abs_array.dat') |
---|
| 129 | !do iT=1,nT |
---|
| 130 | ! do k=1,nS |
---|
| 131 | ! write(113,*), abs_arr(k,iT)*losch**2 |
---|
| 132 | ! end do |
---|
| 133 | !end do |
---|
| 134 | !close(113) |
---|
| 135 | |
---|
[253] | 136 | else |
---|
| 137 | |
---|
[305] | 138 | call bilinearH2H2(wn_arr,temp_arr,abs_arr,wn,temp,abcoef) |
---|
[716] | 139 | abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1 |
---|
[253] | 140 | |
---|
| 141 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
| 142 | ! however our bands are normally thin, so this is no big deal. |
---|
| 143 | |
---|
| 144 | endif |
---|
| 145 | |
---|
| 146 | return |
---|
| 147 | end subroutine interpolateH2H2 |
---|
| 148 | |
---|
| 149 | |
---|
| 150 | !------------------------------------------------------------------------- |
---|
[305] | 151 | subroutine bilinearH2H2(x_arr,y_arr,f2d_arr,x_in,y_in,f) |
---|
[253] | 152 | ! Necessary for interpolation of continuum data |
---|
| 153 | |
---|
| 154 | implicit none |
---|
| 155 | |
---|
| 156 | integer nX,nY,i,j,a,b |
---|
[716] | 157 | parameter(nX=2428) |
---|
| 158 | parameter(nY=10) |
---|
[253] | 159 | |
---|
| 160 | real*8 x_in,y_in,x,y,x1,x2,y1,y2 |
---|
| 161 | real*8 f,f11,f12,f21,f22,fA,fB |
---|
| 162 | real*8 x_arr(nX) |
---|
| 163 | real*8 y_arr(nY) |
---|
| 164 | real*8 f2d_arr(nX,nY) |
---|
| 165 | |
---|
| 166 | integer strlen |
---|
| 167 | character*100 label |
---|
| 168 | label='subroutine bilinear' |
---|
| 169 | |
---|
| 170 | x=x_in |
---|
| 171 | y=y_in |
---|
| 172 | |
---|
[716] | 173 | |
---|
[253] | 174 | ! 1st check we're within the wavenumber range |
---|
| 175 | if ((x.lt.x_arr(2)).or.(x.gt.x_arr(nX-2))) then |
---|
| 176 | f=0.0D+0 |
---|
| 177 | return |
---|
| 178 | else |
---|
| 179 | |
---|
| 180 | ! in the x (wavenumber) direction 1st |
---|
| 181 | i=1 |
---|
| 182 | 10 if (i.lt.(nX+1)) then |
---|
| 183 | if (x_arr(i).gt.x) then |
---|
| 184 | x1=x_arr(i-1) |
---|
| 185 | x2=x_arr(i) |
---|
| 186 | a=i-1 |
---|
| 187 | i=9999 |
---|
| 188 | endif |
---|
| 189 | i=i+1 |
---|
| 190 | goto 10 |
---|
| 191 | endif |
---|
| 192 | endif |
---|
| 193 | |
---|
| 194 | if ((y.lt.y_arr(1)).or.(y.gt.y_arr(nY))) then |
---|
[374] | 195 | write(*,*) 'Warning from bilinearH2H2:' |
---|
[305] | 196 | write(*,*) 'Outside continuum temperature range!' |
---|
[253] | 197 | if(y.lt.y_arr(1))then |
---|
| 198 | y=y_arr(1)+0.01 |
---|
| 199 | endif |
---|
| 200 | if(y.gt.y_arr(nY))then |
---|
| 201 | y=y_arr(nY)-0.01 |
---|
| 202 | endif |
---|
[716] | 203 | else |
---|
[253] | 204 | |
---|
| 205 | ! in the y (temperature) direction 2nd |
---|
| 206 | j=1 |
---|
| 207 | 20 if (j.lt.(nY+1)) then |
---|
| 208 | if (y_arr(j).gt.y) then |
---|
| 209 | y1=y_arr(j-1) |
---|
| 210 | y2=y_arr(j) |
---|
| 211 | b=j-1 |
---|
| 212 | j=9999 |
---|
| 213 | endif |
---|
| 214 | j=j+1 |
---|
| 215 | goto 20 |
---|
| 216 | endif |
---|
[716] | 217 | endif |
---|
[253] | 218 | |
---|
| 219 | f11=f2d_arr(a,b) |
---|
| 220 | f21=f2d_arr(a+1,b) |
---|
| 221 | f12=f2d_arr(a,b+1) |
---|
| 222 | f22=f2d_arr(a+1,b+1) |
---|
[716] | 223 | |
---|
| 224 | call bilinear(f,f11,f21,f12,f22,x,x1,x2,y,y1,y2) |
---|
| 225 | |
---|
[253] | 226 | return |
---|
[305] | 227 | end subroutine bilinearH2H2 |
---|