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

Last change on this file since 2662 was 2662, checked in by dbardet, 3 years ago

Using of is_master keyword to allow nly the master processor to write the messages and avoid too heavy output files

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