source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/interpolateH2H2.F90 @ 1057

Last change on this file since 1057 was 313, checked in by ymipsl, 10 years ago
  • implement splitting of XIOS file for lmdz physics
  • Termination is done properly in parallel by calling MPI_ABORT instead of abort or stop

YM

  • Property svn:executable set to *
File size: 4.2 KB
Line 
1     subroutine interpolateH2H2(wn,temp,pres,abcoef,firstcall,ind)
2
3!==================================================================
4!     
5!     Purpose
6!     -------
7!     Calculates the H2-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 pres               ! pressure               (Pascals)
25
26      ! output
27      double precision abcoef             ! absorption coefficient (m^-1)
28
29      integer nS,nT
30      parameter(nS=2428)
31      parameter(nT=10)
32
33      double precision, parameter :: losch = 2.6867774e19
34      ! Loschmit's number (molecule cm^-3 at STP)
35      ! converts cm^5 molecule^-2 --> cm^-1 amagat^-2
36      ! see Richard et al. 2011, JQSRT for details
37
38      double precision amagat
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
57      integer ind
58
59      if(temp.gt.400)then
60         print*,'Your temperatures are too high for this H2-H2 CIA dataset. If you '
61         print*,'really want to run simulations with hydrogen at T > 400 K, contact'
62         print*,'Robin Wordsworth [rwordsworth@uchicago.edu].'
63         CALL abort_physiq
64      endif
65
66      amagat = (273.15/temp)*(pres/101325.0)
67
68      if(firstcall)then ! called by sugas_corrk only
69         if (is_master) print*,'----------------------------------------------------'
70         if (is_master) print*,'Initialising H2-H2 continuum from HITRAN database...'
71
72!     1.1 Open the ASCII files
73         dt_file=TRIM(datadir)//'/continuum_data/H2-H2_norm_2011.cia'
74
75!$OMP MASTER
76         open(33,file=dt_file,form='formatted',status='old',iostat=ios)
77         if (ios.ne.0) then        ! file not found
78           write(*,*) 'Error from interpolateH2H2'
79           write(*,*) 'Data file ',trim(dt_file),' not found.'
80           write(*,*) 'Check that your path to datagcm:',trim(datadir)
81           write(*,*) 'is correct. You can change it in callphys.def with:'
82           write(*,*) 'datadir = /absolute/path/to/datagcm'
83           write(*,*) 'Also check that the continuum data continuum_data/H2-H2_norm_2011.cia is there.'
84           call abort_physiq
85         else
86
87            do iT=1,nT
88
89               read(33,fmat1) bleh,blah,blah,nres,Ttemp
90               if(nS.ne.nres)then
91                  print*,'Resolution given in file: ',trim(dt_file)
92                  print*,'is ',nres,', which does not match nS.'
93                  print*,'Please adjust nS value in interpolateH2H2.F90'
94                  CALL abort_physiq
95               endif
96               temp_arr(iT)=Ttemp
97
98               do k=1,nS
99                  read(33,*) wn_arr(k),abs_arr(k,it)
100               end do
101
102            end do
103     
104         endif
105         close(33)
106!$OMP END MASTER
107!$OMP BARRIER
108
109         if (is_master) print*,'interpolateH2H2: At wavenumber ',wn,' cm^-1'
110         if (is_master) print*,'   temperature ',temp,' K'
111         if (is_master) print*,'   pressure ',pres,' 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*amagat**2 ! convert to m^-1
121
122         !print*,'We have ',amagat,' amagats of H2'
123         !print*,'So the absorption is ',abcoef,' m^-1'
124
125         ! unlike for Rayleigh scattering, we do not currently weight by the BB function
126         ! however our bands are normally thin, so this is no big deal.
127
128      return
129    end subroutine interpolateH2H2
Note: See TracBrowser for help on using the repository browser.