source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2H2.F90 @ 1157

Last change on this file since 1157 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.

  • Property svn:executable set to *
File size: 4.0 KB
Line 
1     subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall,ind)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2-H2 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
18      implicit none
19
20      ! input
21      double precision wn                 ! wavenumber             (cm^-1)
22      double precision temp               ! temperature            (Kelvin)
23      double precision pres               ! pressure               (Pascals)
24
25      ! output
26      double precision abcoef             ! absorption coefficient (m^-1)
27
28      integer nS,nT
29      parameter(nS=2428)
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      integer ind
57
58      if(temp.gt.400)then
59         print*,'Your temperatures are too high for this H2-H2 CIA dataset. If you '
60         print*,'really want to run simulations with hydrogen at T > 400 K, contact'
61         print*,'Robin Wordsworth [rwordsworth@uchicago.edu].'
62         stop
63      endif
64
65      amagat = (273.15/temp)*(pres/101325.0)
66
67      if(firstcall)then ! called by sugas_corrk only
68         print*,'----------------------------------------------------'
69         print*,'Initialising H2-H2 continuum from HITRAN database...'
70
71!     1.1 Open the ASCII files
72         dt_file=TRIM(datadir)//'/continuum_data/H2-H2_norm_2011.cia'
73
74         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
75         if (ios.ne.0) then        ! file not found
76           write(*,*) 'Error from interpolateH2H2'
77           write(*,*) 'Data file ',trim(dt_file),' not found.'
78           write(*,*) 'Check that your path to datagcm:',trim(datadir)
79           write(*,*) 'is correct. You can change it in callphys.def with:'
80           write(*,*) 'datadir = /absolute/path/to/datagcm'
81           write(*,*) 'Also check that the continuum data continuum_data/H2-H2_norm_2011.cia is there.'
82           call abort
83         else
84
85            do iT=1,nT
86
87               read(33,fmat1) bleh,blah,blah,nres,Ttemp
88               if(nS.ne.nres)then
89                  print*,'Resolution given in file: ',trim(dt_file)
90                  print*,'is ',nres,', which does not match nS.'
91                  print*,'Please adjust nS value in interpolateH2H2.F90'
92                  stop
93               endif
94               temp_arr(iT)=Ttemp
95
96               do k=1,nS
97                  read(33,*) wn_arr(k),abs_arr(k,it)
98               end do
99
100            end do
101     
102         endif
103         close(33)
104
105         print*,'interpolateH2H2: At wavenumber ',wn,' cm^-1'
106         print*,'   temperature ',temp,' K'
107         print*,'   pressure ',pres,' Pa'
108
109      endif
110
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
118         !print*,'We have ',amagat,' amagats of H2'
119         !print*,'So the absorption is ',abcoef,' m^-1'
120
121         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
122         ! however our bands are normally thin, so this is no big deal.
123
124      return
125    end subroutine interpolateH2H2
Note: See TracBrowser for help on using the repository browser.