source: trunk/LMDZ.GENERIC/libf/phystd/interpolateN2H2.F90 @ 1477

Last change on this file since 1477 was 1315, checked in by milmd, 10 years ago

LMDZ.GENERIC. OpenMP directives added in generic physic. When running in pure OpenMP or hybrid OpenMP/MPI, may have some bugs with condense_cloud and wstats routines.

File size: 3.9 KB
Line 
1subroutine interpolateN2H2(wn,temp,presN2,presH2,abcoef,firstcall,ind)
2
3  !==================================================================
4  !     
5  !     Purpose
6  !     -------
7  !     Calculates the N2-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  implicit none
18
19  ! input
20  double precision wn                 ! wavenumber             (cm^-1)
21  double precision temp               ! temperature            (Kelvin)
22  double precision presN2             ! N2 partial pressure    (Pascals)
23  double precision presH2             ! H2 partial pressure    (Pascals)
24
25  ! output
26  double precision abcoef             ! absorption coefficient (m^-1)
27
28  integer nS,nT
29  parameter(nS=1914)
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 amagatN2
38  double precision amagatH2
39  double precision wn_arr(nS)
40  double precision temp_arr(nT)
41  double precision abs_arr(nS,nT)
42
43  integer k,iT
44  logical firstcall
45
46  save wn_arr, temp_arr, abs_arr !read by master
47
48  character*100 dt_file
49  integer strlen,ios
50
51  character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
52
53  character*20 bleh
54  double precision blah, Ttemp
55  integer nres
56  integer ind
57
58  if(temp.gt.400)then
59     print*,'Your temperatures are too high for this N2-H2 CIA dataset.'
60     print*,'Please run mixed N2-H2 atmospheres below T = 400 K.'     
61     stop
62  endif
63
64  amagatN2 = (273.15/temp)*(presN2/101325.0)
65  amagatH2 = (273.15/temp)*(presH2/101325.0)
66
67  if(firstcall)then ! called by sugas_corrk only
68     print*,'----------------------------------------------------'
69     print*,'Initialising N2-H2 continuum from HITRAN database...'
70
71     !     1.1 Open the ASCII files
72     dt_file=TRIM(datadir)//'/continuum_data/N2-H2_2011.cia'
73
74!$OMP MASTER
75     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
76     if (ios.ne.0) then        ! file not found
77        write(*,*) 'Error from interpolateN2H2'
78        write(*,*) 'Data file ',trim(dt_file),' not found.'
79        write(*,*) 'Check that your path to datagcm:',trim(datadir)
80        write(*,*) 'is correct. You can change it in callphys.def with:'
81        write(*,*) 'datadir = /absolute/path/to/datagcm'
82        write(*,*) 'Also check that the continuum data continuum_data/N2-H2_2011.cia is there.'
83        call abort
84     else
85
86        do iT=1,nT
87
88           read(33,fmat1) bleh,blah,blah,nres,Ttemp
89           if(nS.ne.nres)then
90              print*,'Resolution given in file: ',trim(dt_file)
91              print*,'is ',nres,', which does not match nS.'
92              print*,'Please adjust nS value in interpolateN2H2.F90'
93              stop
94           endif
95           temp_arr(iT)=Ttemp
96
97           do k=1,nS
98              read(33,*) wn_arr(k),abs_arr(k,it)
99           end do
100
101        end do
102
103     endif
104     close(33)
105!$OMP END MASTER
106!$OMP BARRIER
107
108     print*,'interpolateN2H2: At wavenumber ',wn,' cm^-1'
109     print*,'   temperature                 ',temp,' K'
110     print*,'   N2 partial pressure         ',presN2,' Pa'
111     print*,'   and H2 partial pressure     ',presH2,' Pa'
112
113  endif
114
115  call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
116
117!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
118!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
119
120  abcoef=abcoef*losch**2*100.0*amagatN2*amagatH2 ! convert to m^-1
121
122!     print*,'We have ',amagatN2,' amagats of N2'
123!     print*,'and     ',amagatH2,' amagats of H2'
124!     print*,'So the absorption is ',abcoef,' m^-1'
125
126
127!  unlike for Rayleigh scattering, we do not currently weight by the BB function
128!  however our bands are normally thin, so this is no big deal.
129
130
131  return
132end subroutine interpolateN2H2
133
Note: See TracBrowser for help on using the repository browser.