[878] | 1 | subroutine interpolateN2N2(wn,temp,pres,abcoef,firstcall,ind) |
---|
[716] | 2 | |
---|
| 3 | !================================================================== |
---|
| 4 | ! |
---|
| 5 | ! Purpose |
---|
| 6 | ! ------- |
---|
| 7 | ! Calculates the N2-N2 CIA opacity, using a lookup table from |
---|
| 8 | ! HITRAN (2011) |
---|
| 9 | ! |
---|
| 10 | ! Authors |
---|
| 11 | ! ------- |
---|
| 12 | ! R. Wordsworth (2011) |
---|
| 13 | ! |
---|
| 14 | !================================================================== |
---|
| 15 | |
---|
| 16 | use datafile_mod, only: datadir |
---|
| 17 | implicit none |
---|
| 18 | |
---|
| 19 | ! input |
---|
| 20 | double precision wn ! wavenumber (cm^-1) |
---|
| 21 | double precision temp ! temperature (Kelvin) |
---|
| 22 | double precision pres ! pressure (Pascals) |
---|
[878] | 23 | integer :: ind |
---|
[716] | 24 | |
---|
| 25 | ! output |
---|
| 26 | double precision abcoef ! absorption coefficient (m^-1) |
---|
| 27 | |
---|
| 28 | integer nS,nT |
---|
| 29 | parameter(nS=582) |
---|
| 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 |
---|
| 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 | |
---|
| 56 | |
---|
| 57 | if(temp.gt.400)then |
---|
| 58 | print*,'Your temperatures are too high for this N2-N2 CIA dataset.' |
---|
| 59 | print*,'Currently, HITRAN provides data for this pair in the range 40-400 K.' |
---|
| 60 | stop |
---|
| 61 | endif |
---|
| 62 | |
---|
| 63 | amagat = (273.15/temp)*(pres/101325.0) |
---|
| 64 | |
---|
| 65 | if(firstcall)then ! called by sugas_corrk only |
---|
| 66 | print*,'----------------------------------------------------' |
---|
| 67 | print*,'Initialising N2-N2 continuum from HITRAN database...' |
---|
| 68 | |
---|
| 69 | ! 1.1 Open the ASCII files |
---|
| 70 | dt_file=TRIM(datadir)//'/continuum_data/N2-N2_2011.cia' |
---|
| 71 | |
---|
| 72 | open(33,file=dt_file,form='formatted',status='old',iostat=ios) |
---|
| 73 | if (ios.ne.0) then ! file not found |
---|
| 74 | write(*,*) 'Error from interpolateN2N2' |
---|
| 75 | write(*,*) 'Data file ',trim(dt_file),' not found.' |
---|
| 76 | write(*,*) 'Check that your path to datagcm:',trim(datadir) |
---|
| 77 | write(*,*) 'is correct. You can change it in callphys.def with:' |
---|
| 78 | write(*,*) 'datadir = /absolute/path/to/datagcm' |
---|
| 79 | write(*,*) 'Also check that the continuum data continuum_data/N2-N2_norm_2011.cia is there.' |
---|
| 80 | call abort |
---|
| 81 | else |
---|
| 82 | |
---|
| 83 | do iT=1,nT |
---|
| 84 | |
---|
| 85 | read(33,fmat1) bleh,blah,blah,nres,Ttemp |
---|
| 86 | if(nS.ne.nres)then |
---|
| 87 | print*,'Resolution given in file: ',trim(dt_file) |
---|
| 88 | print*,'is ',nres,', which does not match nS.' |
---|
| 89 | print*,'Please adjust nS value in interpolateN2N2.F90' |
---|
| 90 | stop |
---|
| 91 | endif |
---|
| 92 | temp_arr(iT)=Ttemp |
---|
| 93 | |
---|
| 94 | do k=1,nS |
---|
| 95 | read(33,*) wn_arr(k),abs_arr(k,it) |
---|
| 96 | end do |
---|
| 97 | |
---|
| 98 | end do |
---|
| 99 | |
---|
| 100 | endif |
---|
| 101 | close(33) |
---|
| 102 | |
---|
| 103 | print*,'interpolateN2N2: At wavenumber ',wn,' cm^-1' |
---|
| 104 | print*,' temperature ',temp,' K' |
---|
| 105 | print*,' pressure ',pres,' Pa' |
---|
| 106 | |
---|
[878] | 107 | endif |
---|
| 108 | call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind) |
---|
[716] | 109 | |
---|
[878] | 110 | ! print*,'the absorption is ',abcoef,' cm^5 molecule^-2' |
---|
| 111 | ! print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2' |
---|
[716] | 112 | |
---|
| 113 | abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1 |
---|
[878] | 114 | ! abcoef=0. |
---|
[716] | 115 | |
---|
[878] | 116 | ! print*,'We have ',amagat,' amagats of N2' |
---|
| 117 | ! print*,'So the absorption is ',abcoef,' m^-1' |
---|
[716] | 118 | |
---|
[878] | 119 | ! unlike for Rayleigh scattering, we do not currently weight by the BB function |
---|
| 120 | ! however our bands are normally thin, so this is no big deal. |
---|
[716] | 121 | |
---|
| 122 | |
---|
| 123 | return |
---|
| 124 | end subroutine interpolateN2N2 |
---|
| 125 | |
---|