source: trunk/LMDZ.GENERIC/libf/phystd/interpolateH2O_self_foreign.F90 @ 2613

Last change on this file since 2613 was 2520, checked in by Guillaume Chaverot, 4 years ago

update of the H2O continuum

File size: 5.8 KB
Line 
1     subroutine interpolateH2O_self_foreign(wn,temp,presS,presF,abcoef,firstcall,ind)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2O continuum self and foreign opacity
8!     using lookup tables from MT_CKD version 3.3
9!
10!     Authors
11!     -------
12!     M. Turbet (2020)
13!     Adapted from the file of R. Wordsworth
14!     
15!==================================================================
16
17      use datafile_mod, only: datadir
18      implicit none
19
20      ! input
21      double precision wn                 ! wavenumber             (cm^-1)
22      double precision temp               ! temperature            (Kelvin)
23      double precision presS              ! self-pressure          (Pascals)
24      double precision presF              ! foreign (air) pressure (Pascals)
25
26      ! output
27      double precision abcoef             ! absorption coefficient (m^-1)
28
29      integer nS,nT, iT
30      parameter(nS=2001)    ! number of wavenumbers
31      parameter(nT=39)      ! number of temperatures
32
33      double precision kB
34      parameter(kB=1.3806488e-23)
35
36      double precision amagatS, amagatF, abcoefS, abcoefF, Nmolec
37      double precision wn_arr(nS)
38      double precision temp_arr(nT)
39      double precision abs_arrS(nS,nT)
40      double precision abs_arrF(nS,nT)
41      double precision data_tmp(nT)
42     
43      character*43 dummy_var1
44      character*46 dummy_var2
45
46      integer nres
47
48      double precision Ttemp
49
50      character(LEN=*), parameter :: fmat1 = "(A43,I6,F6.1)"
51      character(LEN=*), parameter :: fmat2 = "(A46,I6,F6.1)"
52
53
54      integer k,ind
55      logical firstcall
56
57      save wn_arr, temp_arr, abs_arrS, abs_arrF !read by master
58
59      character*100 dt_file
60      integer strlen,ios
61
62      amagatS=(273.15/temp)*(presS/101325.0)
63      amagatF=(273.15/temp)*(presF/101325.0)
64
65      if(firstcall)then ! called by sugas_corrk only
66         print*,'----------------------------------------------------'
67         print*,'Initialising H2O continuum from MT_CKD data...'
68
69!     1.1 Open the ASCII files
70
71!$OMP MASTER
72
73         ! self broadening
74         dt_file=TRIM(datadir)//'/continuum_data/H2O-H2O_continuum_MT_CKD3.3.cia'
75         open(34,file=dt_file,form='formatted',status='old',iostat=ios)
76         if (ios.ne.0) then        ! file not found
77           write(*,*) 'Error from interpolateH2O_cont SELF'
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(*,*) ' Check that there is a H2O-H2O_continuum_MT_CKD3.3.cia'
83           write(*,*)'Continuum file available here:'
84           write(*,*)' https://www.lmd.jussieu.fr/~lmdz/planets/LMDZ.GENERIC/datagcm/continuum_data/'
85           call abort
86         else
87           do iT=1,nT
88             read(34,fmat1) dummy_var1,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*,'Adjust nS value in interpolateH2O_MTCKD...F90'
93               stop
94             endif
95             temp_arr(iT)=Ttemp
96             !write(*,*) 'H2O_H2O'
97             do k=1,nS
98               read(34,*) wn_arr(k),abs_arrS(k,it)
99               !write(*,*) nres,temp_arr(iT),wn_arr(k),abs_arrS(k,it)
100             end do
101           end do
102         endif
103         close(34)
104           
105         ! foreign broadening
106         dt_file=TRIM(datadir)//'/continuum_data/H2O-AIR_continuum_MT_CKD3.3.cia'
107         open(35,file=dt_file,form='formatted',status='old',iostat=ios)
108         if (ios.ne.0) then        ! file not found
109           write(*,*) 'Error from interpolateH2O_cont FOREIGN'
110           write(*,*) 'Data file ',trim(dt_file),' not found.'
111           write(*,*)'Check that your path to datagcm:',trim(datadir)
112           write(*,*)' is correct. You can change it in callphys.def with:'
113           write(*,*)' datadir = /absolute/path/to/datagcm'
114           write(*,*)' Check that there is a H2O-AIR_continuum_MT_CKD3.3.cia'
115           write(*,*)'Continuum file available here:'
116           write(*,*)' https://www.lmd.jussieu.fr/~lmdz/planets/LMDZ.GENERIC/datagcm/continuum_data/'
117           call abort
118         else
119           do iT=1,nT
120             read(35,fmat2) dummy_var2,nres,Ttemp
121             if(nS.ne.nres)then
122               print*,'Resolution given in file: ',trim(dt_file)
123               print*,'is ',nres,', which does not match nS.'
124               print*,'Adjust nS value in interpolateH2O_MTCKD...F90'
125               stop
126             endif
127             temp_arr(iT)=Ttemp
128             !write(*,*) 'H2O_AIR'
129             do k=1,nS
130               read(35,*) wn_arr(k),abs_arrF(k,it)
131               !write(*,*) nres,temp_arr(iT),wn_arr(k),abs_arrF(k,it)
132             end do
133           end do
134         endif
135         close(35)
136
137         print*,'interpolateH2O_MTCKDcont: At wavenumber ',wn,' cm^-1'
138         print*,'   temperature ',temp,' K'
139         print*,'   H2O pressure ',presS,' Pa'
140         print*,'   air pressure ',presF,' Pa'
141!$OMP END MASTER
142!$OMP BARRIER
143         
144      endif
145
146      call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arrS,wn,temp,abcoefS,ind)
147      !print*,'MTCKD - self absorption is ',abcoefS,' cm^2 molecule^-1'
148
149      call bilinearbig(nS,nT,wn_arr,temp_arr,abs_arrF,wn,temp,abcoefF,ind)
150      !print*,'MTCKD - foreign absorption is ',abcoefF,' cm^2 molecule^-1'
151
152! abcoefS and abcoefF are in cm-1 amagat-2
153! First we multiply by amagat^2
154      abcoef = abcoefS*amagatS + abcoefF*amagatF
155      abcoef = abcoef*amagatS
156! Then we convert cm-1 in m-1
157      abcoef = 100.*abcoef
158     
159      !print*,'MTCKD_v3.3 : So the total absorption is ',abcoef,' m^-1'
160      !print*,'for PH2O/PN2/TEMP/wn=',presS,presF,temp,wn
161
162      return
163    end subroutine interpolateH2O_self_foreign
164
165
Note: See TracBrowser for help on using the repository browser.