source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/interpolateH2He.F90 @ 298

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

Less output messages are written. On 20000 cores it is better. In LMDZ, only master of MPI and OpenMP can write.

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