source: trunk/LMDZ.GENERIC/libf/phystd/su_gases.F90 @ 3653

Last change on this file since 3653 was 3653, checked in by gmilcareck, 4 months ago

Generic PCM:
Adding molecular weight for each gas used by the model.
Updating molecular weight for each tracers (to be consistent with gases).
GM

  • Property svn:executable set to *
File size: 5.0 KB
Line 
1subroutine su_gases
2
3  use gases_h
4
5  implicit none
6
7  integer igas, ierr, count
8
9  !==================================================================
10  !     
11  !     Purpose
12  !     -------
13  !     Load atmospheric composition info
14  !     
15  !     Authors
16  !     -------
17  !     R. Wordsworth (2011)
18  !     Allocatable arrays by A. Spiga (2011)
19  !     
20  !==================================================================
21 
22
23!$OMP MASTER
24  ! load gas names from file 'gases.def'
25  open(90,file='gases.def',status='old',form='formatted',iostat=ierr)
26  if (ierr.eq.0) then
27     write(*,*) "sugases.F90: reading file gases.def"
28     read(90,*)
29     read(90,*,iostat=ierr) ngasmx
30     if (ierr.ne.0) then
31        write(*,*) "sugases.F90: error reading number of gases"
32        write(*,*) "   (first line of gases.def) "
33        call abort
34     endif
35
36     print*,ngasmx, " gases found in gases.def. Allocating names and molar fractions..."
37
38     if (.not.allocated(gnom)) allocate(gnom(ngasmx))
39     do igas=1,ngasmx
40        read(90,*,iostat=ierr) gnom(igas)
41        if (ierr.ne.0) then
42           write(*,*) 'sugases.F90: error reading gas names in gases.def...'
43           call abort
44        endif
45     enddo                  !of do igas=1,ngasmx
46
47     vgas=0
48     if(.not.allocated(gfrac)) allocate(gfrac(ngasmx))
49     do igas=1,ngasmx
50        read(90,*,iostat=ierr) gfrac(igas)
51        if (ierr.ne.0) then
52           write(*,*) 'sugases.F90: error reading gas molar fractions in gases.def...'
53           call abort
54        endif
55
56        ! find variable gas (if any)
57        if(gfrac(igas).eq.-1.0)then
58           if(vgas.eq.0)then
59              vgas=igas
60           else
61              print*,'You seem to be choosing two variable gases'
62              print*,'Check that gases.def is correct'
63              call abort
64           endif
65        endif
66
67     enddo                  !of do igas=1,ngasmx
68
69
70     ! assign the 'igas_X' labels
71     count=0
72     IF (.NOT.ALLOCATED(massmol))           ALLOCATE(massmol(ngasmx))
73     do igas=1,ngasmx
74        if (trim(gnom(igas)).eq."H2_" .or. trim(gnom(igas)).eq."H2") then
75           igas_H2=igas
76           massmol(igas_H2)=2.016
77           count=count+1
78        elseif (trim(gnom(igas)).eq."He_" .or. trim(gnom(igas)).eq."He") then
79           igas_He=igas
80           massmol(igas_He)=4.002602
81           count=count+1
82        elseif (trim(gnom(igas)).eq."H2O") then
83           igas_H2O=igas
84           massmol(igas_H2O)=18.01528
85           count=count+1
86        elseif (trim(gnom(igas)).eq."CO2") then
87           igas_CO2=igas
88           massmol(igas_CO2)=44.010
89           count=count+1
90        elseif (trim(gnom(igas)).eq."CO_" .or. trim(gnom(igas)).eq."CO") then
91           igas_CO=igas
92           massmol(igas_CO)=28.010
93           count=count+1
94        elseif (trim(gnom(igas)).eq."N2_" .or. trim(gnom(igas)).eq."N2") then
95           igas_N2=igas
96           massmol(igas_N2)=28.0134
97           count=count+1
98        elseif (trim(gnom(igas)).eq."O2_" .or. trim(gnom(igas)).eq."O2") then
99           igas_O2=igas
100           massmol(igas_O2)=31.998
101           count=count+1
102        elseif (trim(gnom(igas)).eq."SO2") then
103           igas_SO2=igas
104           massmol(igas_SO2)=64.066
105           count=count+1
106        elseif (trim(gnom(igas)).eq."H2S") then
107           igas_H2S=igas
108           massmol(igas_H2S)=34.082
109           count=count+1
110        elseif (trim(gnom(igas)).eq."CH4") then
111           igas_CH4=igas
112           massmol(igas_CH4)=16.0425
113           count=count+1
114        elseif (trim(gnom(igas)).eq."NH3") then
115           igas_NH3=igas
116           massmol(igas_NH3)=17.031
117           count=count+1
118        elseif (trim(gnom(igas)).eq."C2H6") then
119           igas_C2H6=igas
120           massmol(igas_C2H6)=30.069
121           count=count+1
122        elseif (trim(gnom(igas)).eq."C2H2") then
123           igas_C2H2=igas
124           massmol(igas_C2H2)=26.038
125           count=count+1
126        elseif (trim(gnom(igas)).eq."Ar") then
127           igas_Ar=igas
128           massmol(igas_Ar)=39.948
129           count=count+1
130        !! GG MODIF 11 Jan 2019
131       elseif (trim(gnom(igas)).eq."OCS") then
132           igas_OCS=igas
133           massmol(igas_OCS)=60.075
134           count=count+1
135       elseif (trim(gnom(igas)).eq."HCl") then
136           igas_HCl=igas
137           massmol(igas_HCl)=36.461
138           count=count+1
139       elseif (trim(gnom(igas)).eq."HF") then
140           igas_HF=igas
141           massmol(igas_HF)=20.00634
142           count=count+1
143        endif
144     enddo
145
146     if(count.ne.ngasmx)then
147        print*,'Mismatch between ngas and number of recognised gases in sugas_corrk.F90.'
148        print*,'Either we haven`t managed to assign all the gases, or there are duplicates.'
149        print*,'Please try again.'
150     endif
151
152  else
153     write(*,*) 'Cannot find required file "gases.def"'
154     call abort
155  endif
156  close(90)
157!$OMP END MASTER
158!$OMP BARRIER
159
160end subroutine su_gases
Note: See TracBrowser for help on using the repository browser.