[2655] | 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 |
---|
[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 pres ! CH4 partial pressure (Pascals) |
---|
| 26 | |
---|
| 27 | ! output |
---|
| 28 | double precision abcoef ! absorption coefficient (m^-1) |
---|
| 29 | |
---|
| 30 | integer nS,nT |
---|
| 31 | parameter(nS=1018) |
---|
| 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 amagat |
---|
| 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 !read by master |
---|
| 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 |
---|
[2870] | 55 | double precision blah, Ttemp, ztemp |
---|
[2655] | 56 | integer nres |
---|
| 57 | integer ind |
---|
| 58 | |
---|
[2870] | 59 | ztemp = temp |
---|
| 60 | if(ztemp.gt.400)then |
---|
[2655] | 61 | if (strictboundcia) then |
---|
[2662] | 62 | if (is_master) then |
---|
| 63 | print*,'Your temperatures are too high for this CH4-CH4 CIA dataset.' |
---|
| 64 | print*,'Please run mixed CH4-CH4 atmospheres below T = 400 K.' |
---|
| 65 | endif |
---|
[2655] | 66 | stop |
---|
| 67 | else |
---|
[2870] | 68 | !if (is_master) then |
---|
| 69 | ! print*,'Your temperatures are too high for this CH4-CH4 CIA dataset' |
---|
| 70 | ! print*,'you have chosen strictboundcia = ', strictboundcia |
---|
| 71 | ! print*,'*********************************************************' |
---|
| 72 | ! print*,' we allow model to continue but with temp = 400 ' |
---|
| 73 | ! print*,' ... for CH4-CH4 CIA dataset ... ' |
---|
| 74 | ! print*,' ... we assume we know what you are doing ... ' |
---|
| 75 | ! print*,'*********************************************************' |
---|
| 76 | !endif |
---|
| 77 | ztemp = 400 |
---|
[2655] | 78 | endif |
---|
[2870] | 79 | elseif(ztemp.lt.40)then |
---|
[2655] | 80 | if (strictboundcia) then |
---|
[2662] | 81 | if (is_master) then |
---|
| 82 | print*,'Your temperatures are too low for this CH4-CH4 CIA dataset.' |
---|
| 83 | print*,'Please run mixed CH4-CH4 atmospheres above T = 40 K.' |
---|
| 84 | endif |
---|
[2655] | 85 | stop |
---|
| 86 | else |
---|
[2870] | 87 | !if (is_master) then |
---|
| 88 | ! print*,'Your temperatures are too low for this CH4-CH4 CIA dataset' |
---|
| 89 | ! print*,'you have chosen strictboundcia = ', strictboundcia |
---|
| 90 | ! print*,'*********************************************************' |
---|
| 91 | ! print*,' we allow model to continue but with temp = 40 ' |
---|
| 92 | ! print*,' ... for CH4-CH4 CIA dataset ... ' |
---|
| 93 | ! print*,' ... we assume we know what you are doing ... ' |
---|
| 94 | ! print*,'*********************************************************' |
---|
| 95 | !endif |
---|
| 96 | ztemp = 40 |
---|
[2655] | 97 | endif |
---|
| 98 | endif |
---|
| 99 | |
---|
| 100 | amagat = (273.15/temp)*(pres/101325.0) |
---|
| 101 | |
---|
| 102 | if(firstcall)then ! called by sugas_corrk only |
---|
[2662] | 103 | if (is_master) print*,'----------------------------------------------------' |
---|
| 104 | if (is_master) print*,'Initialising CH4-CH4 continuum from HITRAN database...' |
---|
[2655] | 105 | |
---|
| 106 | ! 1.1 Open the ASCII files |
---|
| 107 | dt_file=TRIM(datadir)//'/continuum_data/CH4-CH4_2011.cia' |
---|
| 108 | |
---|
| 109 | !$OMP MASTER |
---|
| 110 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
| 111 | if (ios.ne.0) then ! file not found |
---|
[2662] | 112 | if (is_master) then |
---|
| 113 | write(*,*) 'Error from interpolateCH4CH4' |
---|
| 114 | write(*,*) 'Data file ',trim(dt_file),' not found.' |
---|
| 115 | write(*,*) 'Check that your path to datagcm:',trim(datadir) |
---|
| 116 | write(*,*) 'is correct. You can change it in callphys.def with:' |
---|
| 117 | write(*,*) 'datadir = /absolute/path/to/datagcm' |
---|
| 118 | write(*,*) 'Also check that the continuum data continuum_data/CH4-CH4_2011.cia is there.' |
---|
| 119 | endif |
---|
[2655] | 120 | call abort |
---|
| 121 | else |
---|
| 122 | |
---|
| 123 | do iT=1,nT |
---|
| 124 | |
---|
| 125 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
| 126 | if(nS.ne.nres)then |
---|
[2662] | 127 | if (is_master) then |
---|
| 128 | print*,'Resolution given in file: ',trim(dt_file) |
---|
| 129 | print*,'is ',nres,', which does not match nS.' |
---|
| 130 | print*,'Please adjust nS value in interpolateCH4CH4.F90' |
---|
| 131 | endif |
---|
[2655] | 132 | stop |
---|
| 133 | endif |
---|
| 134 | temp_arr(iT)=Ttemp |
---|
| 135 | |
---|
| 136 | do k=1,nS |
---|
| 137 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
| 138 | end do |
---|
| 139 | |
---|
| 140 | end do |
---|
| 141 | |
---|
| 142 | endif |
---|
| 143 | close(33) |
---|
| 144 | !$OMP END MASTER |
---|
| 145 | !$OMP BARRIER |
---|
[2662] | 146 | if (is_master) then |
---|
| 147 | print*,'interpolateCH4CH4: At wavenumber ',wn,' cm^-1' |
---|
| 148 | print*,' temperature ',temp,' K' |
---|
| 149 | print*,' pressure ',pres,' Pa' |
---|
| 150 | endif |
---|
[2655] | 151 | endif |
---|
| 152 | |
---|
[2870] | 153 | call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,ztemp,abcoef,ind) |
---|
[2655] | 154 | |
---|
| 155 | ! print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
| 156 | ! print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
| 157 | |
---|
| 158 | abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1 |
---|
| 159 | |
---|
| 160 | ! print*,'We have ',amagat,' amagats of CH4' |
---|
| 161 | ! print*,'So the absorption is ',abcoef,' m^-1' |
---|
| 162 | |
---|
| 163 | |
---|
| 164 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
| 165 | ! however our bands are normally thin, so this is no big deal. |
---|
| 166 | |
---|
| 167 | |
---|
| 168 | return |
---|
| 169 | end subroutine interpolateCH4CH4 |
---|
| 170 | |
---|