source: trunk/LMDZ.TITAN/libf/phytitan/su_gases.F90 @ 3529

Last change on this file since 3529 was 1648, checked in by jvatant, 8 years ago

Modifications to custom radiative transfer to Titan
+ Enables an altitude dependant gfrac for CIA computations

-> many radical changes in su_gases and co ..
-> read vertical CH4 profile with call_profilgases
-> Now you need a 'profile.def' that I will add in the deftank

+ Added interpolate CIA routines for CH4
+ Added temporary mean aerosol profile opacity routine (disr_haze)

  • Property svn:executable set to *
File size: 2.5 KB
Line 
1subroutine su_gases(nlayer,tracer)
2
3  use gases_h
4
5  implicit none
6
7  integer,intent(in) :: nlayer
8  logical,intent(in) :: tracer
9 
10  integer igas, ierr
11 
12  !==================================================================
13  !     
14  !     Purpose
15  !     -------
16  !     Load atmospheric composition info
17  !     
18  !     Authors
19  !     -------
20  !     R. Wordsworth (2011)
21  !     Allocatable arrays by A. Spiga (2011)
22  !     Titan's version : J. Vatant d'Ollone (2017)
23  !              * force igas_X labels and def gnom for N2,CH4,H2
24  !              * get rid of variable species ( enrichment dimension will be set in corrk routines )
25  !
26  !==================================================================
27
28!$OMP MASTER
29
30
31  ngasmx   = 3 ! N2, CH4, H2
32
33
34  ! load reference level and reference molar fractions from file 'gases.def'
35  open(90,file='gases.def',status='old',form='formatted',iostat=ierr)
36 
37  if (ierr.eq.0) then
38     write(*,*) "sugases.F90: reading file gases.def"
39     read(90,*) ! header
40
41     ! We allocate gfrac and we set gas molar fractions for the reference level only.
42     ! This will be useful for the cpp_mu and rayleigh routines
43     ! Other gas molar fractions are now set in routine callprofilgases
44     
45     write(*,*) 'sugases.F90: allocating and reading gas molar fractions from reference level in gases.def...'
46     
47     if(.not.allocated(gfrac)) allocate(gfrac(ngasmx,nlayer))
48     
49     read(90,*,iostat=ierr) nivref     
50     if (ierr.ne.0) then
51        write(*,*) "sugases.F90: error reading reference level"
52        write(*,*) "   (first line of gases.def) "
53        call abort
54     endif
55     
56     print*, "layer", nivref, "is reference level found in gases.def ..."
57     
58     do igas=1,ngasmx
59        read(90,*,iostat=ierr) gfrac(igas,nivref)
60        if (ierr.ne.0) then
61           write(*,*) 'sugases.F90: error reading reference gas molar fractions in gases.def... aborting'
62           call abort
63        endif
64     enddo                  !of do igas=1,ngasmx
65
66
67     ! We force gnom = (N2, CH4, H2) and igas_X for Titan
68     
69     igas_N2  = 1
70     igas_CH4 = 2
71     igas_H2  = 3
72
73
74     if (.not.allocated(gnom)) allocate(gnom(ngasmx))
75     gnom(igas_N2)  = "N2_"
76     gnom(igas_CH4) = "CH4"
77     gnom(igas_H2)  = "H2_"
78     
79     
80  else
81     write(*,*) 'Cannot find required file "gases.def"'
82     call abort
83  endif
84 
85  close(90)
86!$OMP END MASTER
87!$OMP BARRIER
88
89end subroutine su_gases
Note: See TracBrowser for help on using the repository browser.