source: trunk/LMDZ.TITAN/libf/phytitan/interpolateN2N2.F90

Last change on this file was 1648, checked in by jvatant, 8 years ago

Modifications to custom radiative transfer to Titan
+ Enables an altitude dependant gfrac for CIA computations

-> many radical changes in su_gases and co ..
-> read vertical CH4 profile with call_profilgases
-> Now you need a 'profile.def' that I will add in the deftank

+ Added interpolate CIA routines for CH4
+ Added temporary mean aerosol profile opacity routine (disr_haze)

File size: 3.7 KB
Line 
1subroutine interpolateN2N2(wn,temp,pres,abcoef,firstcall,ind)
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)
23  integer :: ind
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 !read by master
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!$OMP MASTER
73     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
74     if (ios.ne.0) then        ! file not found
75        write(*,*) 'Error from interpolateN2N2'
76        write(*,*) 'Data file ',trim(dt_file),' not found.'
77        write(*,*) 'Check that your path to datagcm:',trim(datadir)
78        write(*,*) 'is correct. You can change it in callphys.def with:'
79        write(*,*) 'datadir = /absolute/path/to/datagcm'
80        write(*,*) 'Also check that the continuum data continuum_data/N2-N2_2011.cia is there.'
81        call abort
82     else
83
84        do iT=1,nT
85
86           read(33,fmat1) bleh,blah,blah,nres,Ttemp
87           if(nS.ne.nres)then
88              print*,'Resolution given in file: ',trim(dt_file)
89              print*,'is ',nres,', which does not match nS.'
90              print*,'Please adjust nS value in interpolateN2N2.F90'
91              stop
92           endif
93           temp_arr(iT)=Ttemp
94
95           do k=1,nS
96              read(33,*) wn_arr(k),abs_arr(k,it)
97           end do
98
99        end do
100
101     endif
102     close(33)
103!$OMP END MASTER
104!$OMP BARRIER
105
106     print*,'interpolateN2N2: At wavenumber ',wn,' cm^-1'
107     print*,'   temperature ',temp,' K'
108     print*,'   pressure ',pres,' Pa'
109
110  endif
111     call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
112
113!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
114!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
115
116     abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1
117!     abcoef=0.
118
119!     print*,'We have ',amagat,' amagats of N2'
120!     print*,'So the absorption is ',abcoef,' m^-1'
121
122!     unlike for Rayleigh scattering, we do not currently weight by the BB function
123!     however our bands are normally thin, so this is no big deal.
124
125
126  return
127end subroutine interpolateN2N2
128
Note: See TracBrowser for help on using the repository browser.