1 | subroutine inichim_1D(nlayer, nq, pq, qsurf, play, & |
---|
2 | flagh2o,flagthermo) |
---|
3 | |
---|
4 | use tracer_h, only: noms, mmol |
---|
5 | use datafile_mod, only: datadir |
---|
6 | |
---|
7 | implicit none |
---|
8 | |
---|
9 | !======================================================================= |
---|
10 | ! |
---|
11 | ! Purpose: |
---|
12 | ! -------- |
---|
13 | ! |
---|
14 | ! Initialization of the chemistry for use in the building of a new start file |
---|
15 | ! used by program newstart.F |
---|
16 | ! also used by program testphys1d.F |
---|
17 | ! |
---|
18 | ! Authors: |
---|
19 | ! ------- |
---|
20 | ! Initial version 11/2002 by Sebastien Lebonnois |
---|
21 | ! Revised 07/2003 by Monica Angelats-i-Coll to use input files |
---|
22 | ! Modified 10/2008 Identify tracers by their names Ehouarn Millour |
---|
23 | ! Modified 11/2011 Addition of methane Franck Lefevre |
---|
24 | ! Rewritten 04/2012 Franck Lefevre |
---|
25 | ! Rewritten 03/2021 Yassin Jaziri (Use of #Moderntrac-v1 to init thanks traceur.def) |
---|
26 | ! |
---|
27 | ! Arguments: |
---|
28 | ! ---------- |
---|
29 | ! |
---|
30 | ! nlayer Number of atmospheric layers |
---|
31 | ! pq(nlayer,nq) Advected fields, ie chemical species here |
---|
32 | ! qsurf(nq) Amount of tracer on the surface (kg/m2) |
---|
33 | ! play(nlayer) Pressure (Pa) |
---|
34 | ! flagh2o flag for initialisation of h2o (1: yes / 0: no) |
---|
35 | ! flagthermo flag for initialisation of thermosphere only (1: yes / 0: no) |
---|
36 | ! |
---|
37 | !======================================================================= |
---|
38 | |
---|
39 | |
---|
40 | ! inputs : |
---|
41 | |
---|
42 | integer,intent(in) :: nlayer ! Number of atmospheric layers. |
---|
43 | integer,intent(in) :: nq ! number of tracers |
---|
44 | real ,intent(in) :: play(nlayer) ! Mid-layer pressure (Pa). |
---|
45 | integer,intent(in) :: flagh2o ! flag for h2o initialisation |
---|
46 | integer,intent(in) :: flagthermo ! flag for thermosphere initialisation only |
---|
47 | |
---|
48 | ! outputs : |
---|
49 | |
---|
50 | real,intent(out) :: pq(nlayer,nq) ! tracers (kg/kg_of_air) |
---|
51 | real,intent(out) :: qsurf(nq) ! surface values (kg/m2) of tracers |
---|
52 | |
---|
53 | ! local : |
---|
54 | |
---|
55 | real, allocatable :: pf(:) ! pressure in vmr profile files set in traceur.def |
---|
56 | real, allocatable :: qf(:) ! vmr in vmr profile files set in traceur.def |
---|
57 | |
---|
58 | real :: mmean(nlayer) ! mean molecular mass (g) |
---|
59 | real :: pqx(nlayer,nq) ! tracers (vmr) |
---|
60 | real :: qx(nq) ! constant vmr set in traceur.def |
---|
61 | integer :: iq, ilay, iline, nlines, ierr |
---|
62 | |
---|
63 | CHARACTER(len=100) :: qxf(nq) ! vmr profile files set in traceur.def |
---|
64 | CHARACTER(len=100) :: fil ! path files |
---|
65 | character(len=500) :: tracline ! to store a line of text |
---|
66 | |
---|
67 | logical :: foundback = .false. |
---|
68 | |
---|
69 | ! 1. initialization |
---|
70 | |
---|
71 | pq(:,:) = 0. |
---|
72 | qsurf(:) = 0. |
---|
73 | pqx(:,:) = 0. |
---|
74 | qx(:) = 0. |
---|
75 | qxf(:) = 'None' |
---|
76 | nlines = 0 |
---|
77 | |
---|
78 | ! 2. load in traceur.def chemistry data for initialization: |
---|
79 | |
---|
80 | ! Skip nq |
---|
81 | open(90,file='traceur.def',status='old',form='formatted',iostat=ierr) |
---|
82 | if (ierr.eq.0) then |
---|
83 | READ(90,'(A)') tracline |
---|
84 | IF (trim(tracline).eq.'#ModernTrac-v1') THEN ! Test modern traceur.def |
---|
85 | DO |
---|
86 | READ(90,'(A)',iostat=ierr) tracline |
---|
87 | IF (ierr.eq.0) THEN |
---|
88 | IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header |
---|
89 | EXIT |
---|
90 | ENDIF |
---|
91 | ELSE ! If pb, or if reached EOF without having found number of tracer |
---|
92 | write(*,*) "calchim: error reading line of tracers" |
---|
93 | write(*,*) " (first lines of traceur.def) " |
---|
94 | stop |
---|
95 | ENDIF |
---|
96 | ENDDO |
---|
97 | ENDIF ! if modern or standard traceur.def |
---|
98 | else |
---|
99 | write(*,*) "calchim: error opening traceur.def in inichim_1D" |
---|
100 | stop |
---|
101 | endif |
---|
102 | |
---|
103 | ! Get data of tracers |
---|
104 | do iq=1,nq |
---|
105 | read(90,'(A)') tracline |
---|
106 | write(*,*)"inichim_1D: iq=",iq,"noms(iq)=",trim(noms(iq)) |
---|
107 | if (index(tracline,'qx=' ) /= 0) then |
---|
108 | read(tracline(index(tracline,'qx=')+len('qx='):),*) qx(iq) |
---|
109 | write(*,*) ' Parameter value (traceur.def) : qx=', qx(iq) |
---|
110 | else |
---|
111 | write(*,*) ' Parameter value (default) : qx=', qx(iq) |
---|
112 | end if |
---|
113 | if (index(tracline,'qxf=' ) /= 0) then |
---|
114 | read(tracline(index(tracline,'qxf=')+len('qxf='):),*) qxf(iq) |
---|
115 | write(*,*) ' Parameter value (traceur.def) : qxf=', qxf(iq) |
---|
116 | else |
---|
117 | write(*,*) ' Parameter value (default) : qxf=', qxf(iq) |
---|
118 | end if |
---|
119 | end do |
---|
120 | |
---|
121 | close(90) |
---|
122 | |
---|
123 | ! 3. initialization of tracers |
---|
124 | |
---|
125 | ! 3.1 vertical interpolation |
---|
126 | |
---|
127 | do iq=1,nq |
---|
128 | if (qx(iq) /= 0.) then |
---|
129 | pqx(:,iq) = qx(iq) |
---|
130 | else if (qxf(iq) /= 'None') then |
---|
131 | ! Opening file |
---|
132 | fil = trim(datadir)//'/chemical_profiles/'//qxf(iq) |
---|
133 | print*, 'chemical pofile '//trim(noms(iq))//': ', fil |
---|
134 | open(UNIT=90,FILE=fil,STATUS='old',iostat=ierr) |
---|
135 | if (ierr.eq.0) then |
---|
136 | read(90,*) ! read one header line |
---|
137 | do ! get number of lines |
---|
138 | read(90,*,iostat=ierr) |
---|
139 | if (ierr<0) exit |
---|
140 | nlines = nlines + 1 |
---|
141 | end do |
---|
142 | ! allocate reading variable |
---|
143 | allocate(pf(nlines)) |
---|
144 | allocate(qf(nlines)) |
---|
145 | ! read file |
---|
146 | rewind(90) ! restart from the beggining of the file |
---|
147 | read(90,*) ! read one header line |
---|
148 | do iline=1,nlines |
---|
149 | read(90,*) pf(iline), qf(iline) ! pf [Pa], qf [vmr] |
---|
150 | end do |
---|
151 | ! interp in gcm grid |
---|
152 | do ilay=1,nlayer |
---|
153 | call intrplf(log(play(ilay)),pqx(ilay,iq),log(pf),qf,nlines) |
---|
154 | end do |
---|
155 | ! deallocate for next tracer |
---|
156 | deallocate(pf) |
---|
157 | deallocate(qf) |
---|
158 | else |
---|
159 | write(*,*) 'inichim_1D: error opening ', fil |
---|
160 | stop |
---|
161 | endif |
---|
162 | close(90) |
---|
163 | end if |
---|
164 | end do |
---|
165 | |
---|
166 | ! 3.2 background gas |
---|
167 | |
---|
168 | do iq=1,nq |
---|
169 | if (all(pqx(:,iq)==1.)) then |
---|
170 | pqx(:,iq) = 0. |
---|
171 | do ilay=1,nlayer |
---|
172 | pqx(ilay,iq) = 1-sum(pqx(ilay,:)) |
---|
173 | if (pqx(ilay,iq)<=0.) then |
---|
174 | write(*,*) 'inichim_1D: vmr tot > 1 not possible' |
---|
175 | stop |
---|
176 | end if |
---|
177 | end do |
---|
178 | foundback = .true. |
---|
179 | exit ! you found the background gas you can skip others |
---|
180 | end if |
---|
181 | end do |
---|
182 | if (.not.foundback) then |
---|
183 | write(*,*) 'inichim_1D: you need to set a background gas' |
---|
184 | write(*,*) ' by qx=1. in traceur.def' |
---|
185 | stop |
---|
186 | end if |
---|
187 | |
---|
188 | ! 3.3 mmean |
---|
189 | mmean(:) = 0. |
---|
190 | do ilay=1,nlayer |
---|
191 | do iq=1,nq |
---|
192 | mmean(ilay) = mmean(ilay) + pqx(ilay,iq)*mmol(iq) |
---|
193 | end do |
---|
194 | end do |
---|
195 | |
---|
196 | ! 3.4 convert vmr to mmr |
---|
197 | |
---|
198 | do ilay=1,nlayer |
---|
199 | do iq=1,nq |
---|
200 | pq(ilay,iq) = pqx(ilay,iq)*mmol(iq)/mmean(ilay) |
---|
201 | end do |
---|
202 | end do |
---|
203 | |
---|
204 | ! 4. Hard coding |
---|
205 | ! Do whatever you want here to specify pq and qsurf |
---|
206 | ! Or use #ModernTrac-v1 and add another option section 2. |
---|
207 | |
---|
208 | end |
---|