source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2Ocont_CKD.F90 @ 1243

Last change on this file since 1243 was 878, checked in by jleconte, 12 years ago

12/02/2013 == JL

  • Follows previous commit by Aymeric about bilinear interpolations:
    • Extended to all existing continua
    • generalized bilinearbig to work for various size inputs
  • because N2 and H2O continua databases are smaller, improvement around 15% for

an earth case.

File size: 5.9 KB
Line 
1     subroutine interpolateH2Ocont_CKD(wn,temp,presS,presF,abcoef,firstcall,ind)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2O continuum opacity, using a lookup table from
8!     Clough (2005).
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 presS              ! self-pressure          (Pascals)
23      double precision presF              ! foreign (air) pressure (Pascals)
24
25      ! output
26      double precision abcoef             ! absorption coefficient (m^-1)
27
28      integer nS,nT
29      parameter(nS=1001)
30      parameter(nT=11)
31
32      double precision kB
33      parameter(kB=1.3806488e-23)
34
35      double precision amagatS, amagatF, abcoefS, abcoefF, Nmolec
36      double precision wn_arr(nS)
37      double precision temp_arr(nT)
38      double precision abs_arrS(nS,nT)
39      double precision abs_arrF(nS,nT)
40      double precision data_tmp(nT)
41
42      integer k,ind
43      logical firstcall
44
45      save wn_arr, temp_arr, abs_arrS, abs_arrF
46
47      character*100 dt_file
48      integer strlen,ios
49
50      amagatS=(273.15/temp)*(presS/101325.0)
51      amagatF=(273.15/temp)*(presF/101325.0)
52
53      if(firstcall)then ! called by sugas_corrk only
54         print*,'----------------------------------------------------'
55         print*,'Initialising H2O continuum from MT_CKD data...'
56
57!     1.1 Open the ASCII files
58
59         ! nu array
60         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_NU.dat'
61         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
62         if (ios.ne.0) then        ! file not found
63           write(*,*) 'Error from interpolateH2O_cont'
64           write(*,*) 'Data file ',trim(dt_file),' not found.'
65           write(*,*)'Check that your path to datagcm:',trim(datadir)
66           write(*,*)' is correct. You can change it in callphys.def with:'
67           write(*,*)' datadir = /absolute/path/to/datagcm'
68           write(*,*)' Also check that there is a continuum_data/H2O_CONT_NU.dat there.'
69           call abort
70         else
71            do k=1,nS
72               read(33,*) wn_arr(k)
73            enddo
74         endif
75         close(33)
76
77         ! self broadening
78         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_SELF.dat'
79         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
80         if (ios.ne.0) then        ! file not found
81           write(*,*) 'Error from interpolateH2O_cont'
82           write(*,*) 'Data file ',trim(dt_file),' not found.'
83           write(*,*)'Check that your path to datagcm:',trim(datadir)
84           write(*,*)' is correct. You can change it in callphys.def with:'
85           write(*,*)' datadir = /absolute/path/to/datagcm'
86           write(*,*)' Also check that there is a continuum_data/H2O_CONT_SELF.dat there.'
87           call abort
88         else
89            do k=1,nS
90               read(34,*) data_tmp
91               abs_arrS(k,1:nT)=data_tmp(1:nT)
92            end do
93         endif
94         close(34)
95
96         ! foreign (N2+O2+Ar) broadening
97         dt_file=TRIM(datadir)//'/continuum_data/H2O_CONT_FOREIGN.dat'
98         open(35,file=dt_file,form='formatted',status='old',iostat=ios)
99         if (ios.ne.0) then        ! file not found
100           write(*,*) 'Error from interpolateH2O_cont'
101           write(*,*) 'Data file ',trim(dt_file),' not found.'
102           write(*,*)'Check that your path to datagcm:',trim(datadir)
103           write(*,*)' is correct. You can change it in callphys.def with:'
104           write(*,*)' datadir = /absolute/path/to/datagcm'
105           write(*,*)' Also check that there is a continuum_data/H2O_CONT_FOREIGN.dat there.'
106           call abort
107         else
108            do k=1,nS
109               read(35,*) data_tmp
110               abs_arrF(k,1:nT)=data_tmp(1:nT)
111            end do
112         endif
113         close(35)
114
115         temp_arr(1)  = 200.
116         temp_arr(2)  = 250.
117         temp_arr(3)  = 300.
118         temp_arr(4)  = 350.
119         temp_arr(5)  = 400.
120         temp_arr(6)  = 450.
121         temp_arr(7)  = 500.
122         temp_arr(8)  = 550.
123         temp_arr(9)  = 600.
124         temp_arr(10) = 650.
125         temp_arr(11) = 700.
126
127         print*,'interpolateH2Ocont: At wavenumber ',wn,' cm^-1'
128         print*,'   temperature ',temp,' K'
129         print*,'   H2O pressure ',presS,' Pa'
130         print*,'   air pressure ',presF,' Pa'
131         
132      endif
133
134      call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS,ind)
135!      print*,'the self absorption is ',abcoefS,' cm^2 molecule^-1'
136
137      call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF,ind)
138!      print*,'the foreign absorption is ',abcoefF,' cm^2 molecule^-1'
139
140!      print*,'We have ',amagatS,' amagats of H2O vapour'
141!      print*,'and ',amagatF,' amagats of air'
142
143      abcoef = abcoefS*amagatS + abcoefF*amagatF ! Eq. (15) in Clough (1989)
144      abcoef = abcoef*(presS/(presF+presS))      ! take H2O mixing ratio into account
145                                                    ! abs coeffs are given per molecule of H2O
146
147      Nmolec = (presS+presF)/(kB*temp)           ! assume ideal gas
148!      print*,'Total number of molecules per m^3 is',Nmolec
149
150      abcoef = abcoef*Nmolec/(100.0**2)          ! convert to m^-1
151!      print*,'So the total absorption is ',abcoef,' m^-1'
152!      print*,'And optical depth / km : ',1000.0*abcoef
153
154
155      if(wn.gt.500 .and. wn.lt.1400)then
156      elseif(wn.gt.2100 .and. wn.lt.3000)then
157      else
158         abcoef = 0.0
159      endif
160
161      ! unlike for Rayleigh scattering, we do not currently weight by the BB function
162      ! however our bands are normally thin, so this is no big deal.
163
164
165      return
166    end subroutine interpolateH2Ocont_CKD
167
Note: See TracBrowser for help on using the repository browser.