source: trunk/LMDZ.GENERIC/libf/phystd/interpolateCH4CH4.F90 @ 2837

Last change on this file since 2837 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: 5.4 KB
Line 
1subroutine interpolateCH4CH4(wn,temp,pres,abcoef,firstcall,ind)
2
3  !==================================================================
4  !     
5  !     Purpose
6  !     -------
7  !     Calculates the CH4-CH4 CIA opacity, using a lookup table from
8  !     HITRAN (2011)
9  !
10  !     Authors
11  !     -------
12  !     R. Wordsworth (2011)
13  !     
14  !==================================================================
15
16  use callkeys_mod, only: strictboundcia
17  use datafile_mod, only: datadir
18  use mod_phys_lmdz_para, only : is_master
19
20  implicit none
21
22  ! input
23  double precision wn                 ! wavenumber             (cm^-1)
24  double precision temp               ! temperature            (Kelvin)
25  double precision pres               ! CH4 partial pressure    (Pascals)
26
27  ! output
28  double precision abcoef             ! absorption coefficient (m^-1)
29
30  integer nS,nT
31  parameter(nS=1018)
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 amagat
40  double precision wn_arr(nS)
41  double precision temp_arr(nT)
42  double precision abs_arr(nS,nT)
43
44  integer k,iT
45  logical firstcall
46
47  save wn_arr, temp_arr, abs_arr !read by master
48
49  character*100 dt_file
50  integer strlen,ios
51
52  character(LEN=*), parameter :: fmat1 = "(A20,F10.3,F10.3,I7,F7.1,E10.3,F5.3)"
53
54  character*20 bleh
55  double precision blah, Ttemp
56  integer nres
57  integer ind
58
59  if(temp.gt.400)then
60    if (strictboundcia) then
61      if (is_master) then
62        print*,'Your temperatures are too high for this CH4-CH4 CIA dataset.'
63        print*,'Please run mixed CH4-CH4 atmospheres below T = 400 K.'     
64      endif
65      stop
66    else
67      if (is_master) then
68        print*,'Your temperatures are too high for this CH4-CH4 CIA dataset'
69        print*,'you have chosen strictboundcia = ', strictboundcia
70        print*,'*********************************************************'
71        print*,' we allow model to continue but with temp = 400          '
72        print*,'  ...       for CH4-CH4 CIA dataset        ...           '
73        print*,'  ... we assume we know what you are doing ...           '
74        print*,'*********************************************************'
75      endif
76      temp = 400
77    endif
78  elseif(temp.lt.40)then
79    if (strictboundcia) then
80      if (is_master) then
81        print*,'Your temperatures are too low for this CH4-CH4 CIA dataset.'
82        print*,'Please run mixed CH4-CH4 atmospheres above T = 40 K.'     
83      endif
84      stop
85    else
86      if (is_master) then
87        print*,'Your temperatures are too low for this CH4-CH4 CIA dataset'
88        print*,'you have chosen strictboundcia = ', strictboundcia
89        print*,'*********************************************************'
90        print*,' we allow model to continue but with temp = 40           '
91        print*,'  ...       for CH4-CH4 CIA dataset        ...           '
92        print*,'  ... we assume we know what you are doing ...           '
93        print*,'*********************************************************'
94      endif
95      temp = 40
96    endif
97  endif
98
99  amagat = (273.15/temp)*(pres/101325.0)
100
101  if(firstcall)then ! called by sugas_corrk only
102     if (is_master) print*,'----------------------------------------------------'
103     if (is_master) print*,'Initialising CH4-CH4 continuum from HITRAN database...'
104
105     !     1.1 Open the ASCII files
106     dt_file=TRIM(datadir)//'/continuum_data/CH4-CH4_2011.cia'
107
108!$OMP MASTER
109     open(33,file=dt_file,form='formatted',status='old',iostat=ios)
110     if (ios.ne.0) then        ! file not found
111        if (is_master) then
112          write(*,*) 'Error from interpolateCH4CH4'
113          write(*,*) 'Data file ',trim(dt_file),' not found.'
114          write(*,*) 'Check that your path to datagcm:',trim(datadir)
115          write(*,*) 'is correct. You can change it in callphys.def with:'
116          write(*,*) 'datadir = /absolute/path/to/datagcm'
117          write(*,*) 'Also check that the continuum data continuum_data/CH4-CH4_2011.cia is there.'
118        endif
119        call abort
120     else
121
122        do iT=1,nT
123
124           read(33,fmat1) bleh,blah,blah,nres,Ttemp
125           if(nS.ne.nres)then
126              if (is_master) then
127                print*,'Resolution given in file: ',trim(dt_file)
128                print*,'is ',nres,', which does not match nS.'
129                print*,'Please adjust nS value in interpolateCH4CH4.F90'
130              endif
131              stop
132           endif
133           temp_arr(iT)=Ttemp
134
135           do k=1,nS
136              read(33,*) wn_arr(k),abs_arr(k,it)
137           end do
138
139        end do
140
141     endif
142     close(33)
143!$OMP END MASTER
144!$OMP BARRIER
145     if (is_master) then
146       print*,'interpolateCH4CH4: At wavenumber ',wn,' cm^-1'
147       print*,'   temperature                 ',temp,' K'
148       print*,'   pressure                    ',pres,' Pa'
149     endif
150  endif
151
152  call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arr,wn,temp,abcoef,ind)
153
154!     print*,'the absorption is ',abcoef,' cm^5 molecule^-2'
155!     print*,'or ',abcoef*losch**2,' cm^-1 amagat^-2'
156
157  abcoef=abcoef*losch**2*100.0*amagat**2 ! convert to m^-1
158
159!     print*,'We have ',amagat,' amagats of CH4'
160!     print*,'So the absorption is ',abcoef,' m^-1'
161
162
163!  unlike for Rayleigh scattering, we do not currently weight by the BB function
164!  however our bands are normally thin, so this is no big deal.
165
166
167  return
168end subroutine interpolateCH4CH4
169
Note: See TracBrowser for help on using the repository browser.