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