[878] | 1 | subroutine interpolateH2He(wn,temp,presH2,presHe,abcoef,firstcall,ind) |
---|
[716] | 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 | |
---|
[2655] | 16 | use callkeys_mod, only: H2orthopara_mixture, strictboundcia |
---|
[716] | 17 | use datafile_mod, only: datadir |
---|
[2662] | 18 | use mod_phys_lmdz_para, only : is_master |
---|
[873] | 19 | |
---|
[716] | 20 | implicit none |
---|
| 21 | |
---|
| 22 | ! input |
---|
| 23 | double precision wn ! wavenumber (cm^-1) |
---|
| 24 | double precision temp ! temperature (Kelvin) |
---|
| 25 | double precision presH2 ! H2 partial pressure (Pascals) |
---|
| 26 | double precision presHe ! He partial pressure (Pascals) |
---|
| 27 | |
---|
| 28 | ! output |
---|
| 29 | double precision abcoef ! absorption coefficient (m^-1) |
---|
| 30 | |
---|
| 31 | integer nS,nT |
---|
| 32 | parameter(nT=10) |
---|
| 33 | |
---|
| 34 | double precision, parameter :: losch = 2.6867774e19 |
---|
| 35 | ! Loschmit's number (molecule cm^-3 at STP) |
---|
| 36 | ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2 |
---|
| 37 | ! see Richard et al. 2011, JQSRT for details |
---|
| 38 | |
---|
| 39 | double precision amagatH2 |
---|
| 40 | double precision amagatHe |
---|
| 41 | double precision temp_arr(nT) |
---|
[2667] | 42 | double precision, dimension(:), allocatable :: wn_arr |
---|
| 43 | double precision, dimension(:,:), allocatable :: abs_arr |
---|
[716] | 44 | |
---|
| 45 | integer k,iT |
---|
| 46 | logical firstcall |
---|
| 47 | |
---|
[1315] | 48 | save wn_arr, temp_arr, abs_arr !read by master |
---|
[716] | 49 | |
---|
| 50 | character*100 dt_file |
---|
| 51 | integer strlen,ios |
---|
| 52 | |
---|
| 53 | character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)" |
---|
| 54 | |
---|
| 55 | character*20 bleh |
---|
| 56 | double precision blah, Ttemp |
---|
| 57 | integer nres |
---|
| 58 | |
---|
[878] | 59 | integer ind |
---|
[716] | 60 | |
---|
[2667] | 61 | |
---|
| 62 | if (H2orthopara_mixture .eq. "hot") then ! .and. (temp .gt. 3000.0)) then |
---|
| 63 | ! print*,"We're in the Hot Jupiter case " |
---|
| 64 | nS = 19981 |
---|
| 65 | if (temp .gt. 9900) then |
---|
| 66 | if (strictboundcia) then |
---|
| 67 | if (is_master) then |
---|
| 68 | print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case).' |
---|
| 69 | print*,'Please run mixed H2-He atmospheres below T = 9900.0 K.' |
---|
| 70 | endif !is_master |
---|
| 71 | stop |
---|
| 72 | else |
---|
| 73 | if (is_master) then |
---|
| 74 | print*,'Your temperatures are too high for this H2-He CIA dataset (Hot Jupiter case)' |
---|
| 75 | print*,'you have chosen strictboundcia = ', strictboundcia |
---|
| 76 | print*,'*********************************************************' |
---|
| 77 | print*,' we allow model to continue but with temp = 9900.0 K ' |
---|
| 78 | print*,' ... for H2-He CIA dataset ... ' |
---|
| 79 | print*,' ... we assume we know what you are doing ... ' |
---|
| 80 | print*,'*********************************************************' |
---|
| 81 | endif !is_master |
---|
| 82 | temp = 9900.0 |
---|
| 83 | endif ! of stricbound cia |
---|
| 84 | endif ! temp .gt. 9900 |
---|
| 85 | else !not in Hot Jupiter case |
---|
| 86 | nS = 2428 |
---|
| 87 | if(temp.gt.400)then |
---|
| 88 | if (strictboundcia) then |
---|
| 89 | if (is_master) then |
---|
| 90 | print*,'Your temperatures are too high for this H2-He CIA dataset.' |
---|
| 91 | print*,'Please run mixed H2-He atmospheres below T = 400 K.' |
---|
| 92 | endif ! is_master |
---|
| 93 | stop |
---|
| 94 | else |
---|
| 95 | if (is_master) then |
---|
| 96 | print*,'Your temperatures are too high for this H2-He CIA dataset' |
---|
| 97 | print*,'you have chosen strictboundcia = ', strictboundcia |
---|
| 98 | print*,'*********************************************************' |
---|
| 99 | print*,' we allow model to continue but with temp = 400 ' |
---|
| 100 | print*,' ... for H2-He CIA dataset ... ' |
---|
| 101 | print*,' ... we assume we know what you are doing ... ' |
---|
| 102 | print*,'*********************************************************' |
---|
| 103 | endif ! is_master |
---|
| 104 | temp = 400 |
---|
| 105 | endif !of strictboundcia |
---|
| 106 | elseif(temp.lt.40)then |
---|
| 107 | if (strictboundcia) then |
---|
| 108 | if (is_master) then |
---|
| 109 | print*,'Your temperatures are too low for this H2-He CIA dataset.' |
---|
| 110 | print*,'Please run mixed H2-He atmospheres above T = 40 K.' |
---|
| 111 | endif ! is_master |
---|
| 112 | stop |
---|
| 113 | else |
---|
| 114 | if (is_master) then |
---|
| 115 | print*,'Your temperatures are too low for this H2-He CIA dataset' |
---|
| 116 | print*,'you have chosen strictboundcia = ', strictboundcia |
---|
| 117 | print*,'*********************************************************' |
---|
| 118 | print*,' we allow model to continue but with temp = 40 ' |
---|
| 119 | print*,' ... for H2-He CIA dataset ... ' |
---|
| 120 | print*,' ... we assume we know what you are doing ... ' |
---|
| 121 | print*,'*********************************************************' |
---|
| 122 | endif ! is_master |
---|
| 123 | temp = 40 |
---|
| 124 | endif !of strictboundcia |
---|
| 125 | endif ! of (temp .gt. 400) |
---|
| 126 | endif ! of (H2orthopara_mixture .eq. "hot") |
---|
| 127 | |
---|
| 128 | if (.not. allocated(wn_arr)) allocate(wn_arr(nS)) |
---|
| 129 | if (.not. allocated(abs_arr)) allocate(abs_arr(nS,nT)) |
---|
| 130 | |
---|
[716] | 131 | amagatH2 = (273.15/temp)*(presH2/101325.0) |
---|
| 132 | amagatHe = (273.15/temp)*(presHe/101325.0) |
---|
| 133 | |
---|
| 134 | if(firstcall)then ! called by sugas_corrk only |
---|
[2662] | 135 | if (is_master) print*,'----------------------------------------------------' |
---|
| 136 | if (is_master) print*,'Initialising H2-He continuum from HITRAN database...' |
---|
[716] | 137 | |
---|
| 138 | ! 1.1 Open the ASCII files |
---|
[2655] | 139 | if (H2orthopara_mixture.eq."normal") then |
---|
| 140 | dt_file=TRIM(datadir)//'/continuum_data/H2-He_norm_2011.cia' |
---|
| 141 | else if (H2orthopara_mixture.eq."equilibrium") then |
---|
| 142 | dt_file=TRIM(datadir)//'/continuum_data/H2-He_eq_2011.cia' |
---|
[2667] | 143 | else if (H2orthopara_mixture .eq. "hot") then !we use a dataset than can go as high as 9900 K |
---|
| 144 | dt_file=TRIM(datadir)//'/continuum_data/H2-He_2011.cia' |
---|
[2655] | 145 | endif |
---|
[1315] | 146 | |
---|
| 147 | !$OMP MASTER |
---|
[716] | 148 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
| 149 | if (ios.ne.0) then ! file not found |
---|
[2662] | 150 | if (is_master) then |
---|
| 151 | write(*,*) 'Error from interpolateH2He' |
---|
| 152 | write(*,*) 'Data file ',trim(dt_file),' not found.' |
---|
| 153 | write(*,*) 'Check that your path to datagcm:',trim(datadir) |
---|
| 154 | write(*,*) 'is correct. You can change it in callphys.def with:' |
---|
| 155 | write(*,*) 'datadir = /absolute/path/to/datagcm' |
---|
| 156 | write(*,*) 'Also check that the continuum data is there.' |
---|
| 157 | endif |
---|
[716] | 158 | call abort |
---|
| 159 | else |
---|
| 160 | |
---|
| 161 | do iT=1,nT |
---|
| 162 | |
---|
| 163 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
| 164 | if(nS.ne.nres)then |
---|
[2662] | 165 | if (is_master) then |
---|
| 166 | print*,'Resolution given in file: ',trim(dt_file) |
---|
| 167 | print*,'is ',nres,', which does not match nS.' |
---|
| 168 | print*,'Please adjust nS value in interpolateH2He.F90' |
---|
| 169 | endif |
---|
[716] | 170 | stop |
---|
| 171 | endif |
---|
| 172 | temp_arr(iT)=Ttemp |
---|
| 173 | |
---|
| 174 | do k=1,nS |
---|
| 175 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
| 176 | end do |
---|
| 177 | |
---|
| 178 | end do |
---|
| 179 | |
---|
| 180 | endif |
---|
| 181 | close(33) |
---|
[1315] | 182 | !$OMP END MASTER |
---|
| 183 | !$OMP BARRIER |
---|
[2662] | 184 | if (is_master) then |
---|
| 185 | print*,'interpolateH2He: At wavenumber ',wn,' cm^-1' |
---|
| 186 | print*,' temperature ',temp,' K' |
---|
| 187 | print*,' H2 partial pressure ',presH2,' Pa' |
---|
| 188 | print*,' and He partial pressure ',presHe,' Pa' |
---|
| 189 | endif |
---|
[873] | 190 | endif |
---|
[716] | 191 | |
---|
[878] | 192 | call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind) |
---|
[716] | 193 | |
---|
[873] | 194 | !print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
| 195 | !print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
| 196 | |
---|
[716] | 197 | abcoef=abcoef*losch**2*100.0*amagatH2*amagatHe ! convert to m^-1 |
---|
| 198 | |
---|
[873] | 199 | !print*,'We have ',amagatH2,' amagats of H2' |
---|
| 200 | !print*,'and ',amagatHe,' amagats of He' |
---|
| 201 | !print*,'So the absorption is ',abcoef,' m^-1' |
---|
[716] | 202 | |
---|
| 203 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
| 204 | ! however our bands are normally thin, so this is no big deal. |
---|
| 205 | |
---|
| 206 | return |
---|
| 207 | end subroutine interpolateH2He |
---|