source: trunk/LMDZ.GENERIC/libf/phystd/initracer.F @ 2326

Last change on this file since 2326 was 2278, checked in by emillour, 5 years ago

Generic GCM:

  • Get rid of "zerophys"; modern Fortran syntax is better to initialize an array.

EM

File size: 11.2 KB
Line 
1      SUBROUTINE initracer(ngrid,nq,nametrac)
2
3      use surfdat_h, ONLY: dryness, watercaptag
4      USE tracer_h
5      USE callkeys_mod, only: water
6      IMPLICIT NONE
7c=======================================================================
8c   subject:
9c   --------
10c   Initialization related to tracer
11c   (transported dust, water, chemical species, ice...)
12c
13c   Name of the tracer
14c
15c   Test of dimension :
16c   Initialize COMMON tracer in tracer.h, using tracer names provided
17c   by the argument nametrac
18c
19c   author: F.Forget
20c   ------
21c            Ehouarn Millour (oct. 2008) identify tracers by their names
22c=======================================================================
23
24      integer,intent(in) :: ngrid,nq
25      character(len=30),intent(in) :: nametrac(nq) ! name of the tracer from dynamics
26
27      character(len=30) :: txt ! to store some text
28      integer iq,ig,count
29      real r0_lift , reff_lift
30
31
32
33c-----------------------------------------------------------------------
34c  radius(nq)      ! aerosol particle radius (m)
35c  rho_q(nq)       ! tracer densities (kg.m-3)
36c  qext(nq)        ! Single Scat. Extinction coeff at 0.67 um
37c  alpha_lift(nq)  ! saltation vertical flux/horiz flux ratio (m-1)
38c  alpha_devil(nq) ! lifting coeeficient by dust devil
39c  rho_dust          ! Mars dust density
40c  rho_ice           ! Water ice density
41c  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
42c  varian            ! Characteristic variance of log-normal distribution
43c-----------------------------------------------------------------------
44
45       nqtot=nq
46       !! we allocate once for all arrays in common in tracer_h.F90
47       !! (supposedly those are not used before call to initracer)
48       IF (.NOT.ALLOCATED(noms))         ALLOCATE(noms(nq))
49       IF (.NOT.ALLOCATED(mmol))         ALLOCATE(mmol(nq))
50       IF (.NOT.ALLOCATED(radius))       ALLOCATE(radius(nq))
51       IF (.NOT.ALLOCATED(rho_q))        ALLOCATE(rho_q(nq))
52       IF (.NOT.ALLOCATED(qext))         ALLOCATE(qext(nq))
53       IF (.NOT.ALLOCATED(alpha_lift))   ALLOCATE(alpha_lift(nq))
54       IF (.NOT.ALLOCATED(alpha_devil))  ALLOCATE(alpha_devil(nq))
55       IF (.NOT.ALLOCATED(qextrhor))     ALLOCATE(qextrhor(nq))
56       IF (.NOT.ALLOCATED(igcm_dustbin)) ALLOCATE(igcm_dustbin(nq))
57       !! initialization
58       alpha_lift(:)=0.
59       alpha_devil(:)=0.
60       
61       ! Added by JVO 2017 : these arrays are handled later
62       ! -> initialization is the least we can do, please !!!
63       radius(:)=0.
64       qext(:)=0.
65
66
67! Initialization: copy tracer names from dynamics
68        do iq=1,nq
69          noms(iq)=nametrac(iq)
70          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
71        enddo
72
73
74! Identify tracers by their names: (and set corresponding values of mmol)
75      ! 0. initialize tracer indexes to zero:
76      ! NB: igcm_* indexes are commons in 'tracer.h'
77      do iq=1,nq
78        igcm_dustbin(iq)=0
79      enddo
80      igcm_dust_mass=0
81      igcm_dust_number=0
82      igcm_h2o_vap=0
83      igcm_h2o_ice=0
84      igcm_co2=0
85      igcm_co=0
86      igcm_o=0
87      igcm_o1d=0
88      igcm_o2=0
89      igcm_o3=0
90      igcm_h=0
91      igcm_h2=0
92      igcm_oh=0
93      igcm_ho2=0
94      igcm_h2o2=0
95      igcm_n2=0
96      igcm_n=0
97      igcm_n2d=0
98      igcm_no=0
99      igcm_no2=0
100      igcm_ar=0
101      igcm_ar_n2=0
102      igcm_co2_ice=0
103
104      igcm_ch4=0
105      igcm_ch3=0
106      igcm_ch=0
107      igcm_3ch2=0
108      igcm_1ch2=0
109      igcm_cho=0
110      igcm_ch2o=0
111      igcm_ch3o=0
112      igcm_c=0
113      igcm_c2=0
114      igcm_c2h=0
115      igcm_c2h2=0
116      igcm_c2h3=0
117      igcm_c2h4=0
118      igcm_c2h6=0
119      igcm_ch2co=0
120      igcm_ch3co=0
121      igcm_hcaer=0
122
123
124
125      ! 1. find dust tracers
126      count=0
127
128      ! 2. find chemistry and water tracers
129      do iq=1,nq
130        if (noms(iq).eq."co2") then
131          igcm_co2=iq
132          mmol(igcm_co2)=44.
133          count=count+1
134!          write(*,*) 'co2: count=',count
135        endif
136        if (noms(iq).eq."co2_ice") then
137          igcm_co2_ice=iq
138          mmol(igcm_co2_ice)=44.
139          count=count+1
140!          write(*,*) 'co2_ice: count=',count
141        endif
142        if (noms(iq).eq."h2o_vap") then
143          igcm_h2o_vap=iq
144          mmol(igcm_h2o_vap)=18.
145          count=count+1
146!          write(*,*) 'h2o_vap: count=',count
147        endif
148        if (noms(iq).eq."h2o_ice") then
149          igcm_h2o_ice=iq
150          mmol(igcm_h2o_ice)=18.
151          count=count+1
152!          write(*,*) 'h2o_ice: count=',count
153        endif
154        if (noms(iq).eq."co") then
155          igcm_co=iq
156          mmol(igcm_co)=28.
157          count=count+1
158        endif
159        if (noms(iq).eq."o") then
160          igcm_o=iq
161          mmol(igcm_o)=16.
162          count=count+1
163        endif
164        if (noms(iq).eq."o1d") then
165          igcm_o1d=iq
166          mmol(igcm_o1d)=16.
167          count=count+1
168        endif
169        if (noms(iq).eq."o2") then
170          igcm_o2=iq
171          mmol(igcm_o2)=32.
172          count=count+1
173        endif
174        if (noms(iq).eq."o3") then
175          igcm_o3=iq
176          mmol(igcm_o3)=48.
177          count=count+1
178        endif
179        if (noms(iq).eq."h") then
180          igcm_h=iq
181          mmol(igcm_h)=1.
182          count=count+1
183        endif
184        if (noms(iq).eq."h2") then
185          igcm_h2=iq
186          mmol(igcm_h2)=2.
187          count=count+1
188        endif
189        if (noms(iq).eq."oh") then
190          igcm_oh=iq
191          mmol(igcm_oh)=17.
192          count=count+1
193        endif
194        if (noms(iq).eq."ho2") then
195          igcm_ho2=iq
196          mmol(igcm_ho2)=33.
197          count=count+1
198        endif
199        if (noms(iq).eq."h2o2") then
200          igcm_h2o2=iq
201          mmol(igcm_h2o2)=34.
202          count=count+1
203        endif
204        if (noms(iq).eq."n2") then
205          igcm_n2=iq
206          mmol(igcm_n2)=28.
207          count=count+1
208        endif
209        if (noms(iq).eq."ch4") then
210          igcm_ch4=iq
211          mmol(igcm_ch4)=16.
212          count=count+1
213        endif
214        if (noms(iq).eq."ar") then
215          igcm_ar=iq
216          mmol(igcm_ar)=40.
217          count=count+1
218        endif
219        if (noms(iq).eq."n") then
220          igcm_n=iq
221          mmol(igcm_n)=14.
222          count=count+1
223        endif
224        if (noms(iq).eq."no") then
225          igcm_no=iq
226          mmol(igcm_no)=30.
227          count=count+1
228        endif
229        if (noms(iq).eq."no2") then
230          igcm_no2=iq
231          mmol(igcm_no2)=46.
232          count=count+1
233        endif
234        if (noms(iq).eq."n2d") then
235          igcm_n2d=iq
236          mmol(igcm_n2d)=28.
237          count=count+1
238        endif
239        if (noms(iq).eq."ch3") then
240          igcm_ch3=iq
241          mmol(igcm_ch3)=15.
242          count=count+1
243        endif
244        if (noms(iq).eq."ch") then
245          igcm_ch=iq
246          mmol(igcm_ch)=13.
247          count=count+1
248        endif
249        if (noms(iq).eq."3ch2") then
250          igcm_3ch2=iq
251          mmol(igcm_3ch2)=14.
252          count=count+1
253        endif
254        if (noms(iq).eq."1ch2") then
255          igcm_1ch2=iq
256          mmol(igcm_1ch2)=14.
257          count=count+1
258        endif
259        if (noms(iq).eq."cho") then
260          igcm_cho=iq
261          mmol(igcm_cho)=29.
262          count=count+1
263        endif
264        if (noms(iq).eq."ch2o") then
265          igcm_ch2o=iq
266          mmol(igcm_ch2o)=30.
267          count=count+1
268        endif
269        if (noms(iq).eq."ch3o") then
270          igcm_ch3o=iq
271          mmol(igcm_ch3o)=31.
272          count=count+1
273        endif
274        if (noms(iq).eq."c") then
275          igcm_c=iq
276          mmol(igcm_c)=12.
277          count=count+1
278        endif
279        if (noms(iq).eq."c2") then
280          igcm_c2=iq
281          mmol(igcm_c2)=24.
282          count=count+1
283        endif
284        if (noms(iq).eq."c2h") then
285          igcm_c2h=iq
286          mmol(igcm_c2h)=25.
287          count=count+1
288        endif
289        if (noms(iq).eq."c2h2") then
290          igcm_c2h2=iq
291          mmol(igcm_c2h2)=26.
292          count=count+1
293        endif
294        if (noms(iq).eq."c2h3") then
295          igcm_c2h3=iq
296          mmol(igcm_c2h3)=27.
297          count=count+1
298        endif
299        if (noms(iq).eq."c2h4") then
300          igcm_c2h4=iq
301          mmol(igcm_c2h4)=28.
302          count=count+1
303        endif
304        if (noms(iq).eq."c2h6") then
305          igcm_c2h6=iq
306          mmol(igcm_c2h6)=30.
307          count=count+1
308        endif
309        if (noms(iq).eq."ch2co") then
310          igcm_ch2co=iq
311          mmol(igcm_ch2co)=42.
312          count=count+1
313        endif
314        if (noms(iq).eq."ch3co") then
315          igcm_ch3co=iq
316          mmol(igcm_ch3co)=43.
317          count=count+1
318        endif
319        if (noms(iq).eq."hcaer") then
320          igcm_hcaer=iq
321          mmol(igcm_hcaer)=50.
322          count=count+1
323        endif
324      enddo ! of do iq=1,nq
325     
326      ! check that we identified all tracers:
327      if (count.ne.nq) then
328        write(*,*) "initracer: found only ",count," tracers"
329        write(*,*) "               expected ",nq
330        do iq=1,count
331          write(*,*)'      ',iq,' ',trim(noms(iq))
332        enddo
333!        stop
334      else
335        write(*,*) "initracer: found all expected tracers, namely:"
336        do iq=1,nq
337          write(*,*)'      ',iq,' ',trim(noms(iq))
338        enddo
339      endif
340
341
342c------------------------------------------------------------
343c     Initialisation tracers ....
344c------------------------------------------------------------
345      rho_q(1:nq)=0
346
347      rho_dust=2500.  ! Mars dust density (kg.m-3)
348      rho_ice=920.    ! Water ice density (kg.m-3)
349      rho_co2=1620.   ! CO2 ice density (kg.m-3)
350
351c     Initialization for water vapor
352c     ------------------------------
353      if(water) then
354         radius(igcm_h2o_vap)=0.
355         Qext(igcm_h2o_vap)=0.
356         alpha_lift(igcm_h2o_vap) =0.
357         alpha_devil(igcm_h2o_vap)=0.
358         qextrhor(igcm_h2o_vap)= 0.
359
360         !! this is defined in surfdat_h.F90
361         IF (.not.ALLOCATED(dryness)) ALLOCATE(dryness(ngrid))
362         IF (.not.ALLOCATED(watercaptag)) ALLOCATE(watercaptag(ngrid))
363
364         do ig=1,ngrid
365           if (ngrid.ne.1) watercaptag(ig)=.false.
366           dryness(ig) = 1.
367         enddo
368
369
370           radius(igcm_h2o_ice)=3.e-6
371           rho_q(igcm_h2o_ice)=rho_ice
372           Qext(igcm_h2o_ice)=0.
373!           alpha_lift(igcm_h2o_ice) =0.
374!           alpha_devil(igcm_h2o_ice)=0.
375           qextrhor(igcm_h2o_ice)= (3./4.)*Qext(igcm_h2o_ice)
376     $       / (rho_ice*radius(igcm_h2o_ice))
377
378
379      end if  ! (water)
380
381
382!
383!     some extra (possibly redundant) sanity checks for tracers:
384!     ---------------------------------------------------------
385       if (water) then
386       ! verify that we indeed have h2o_vap and h2o_ice tracers
387         if (igcm_h2o_vap.eq.0) then
388           write(*,*) "initracer: error !!"
389           write(*,*) "  cannot use water option without ",
390     &                "an h2o_vap tracer !"
391           stop
392         endif
393         if (igcm_h2o_ice.eq.0) then
394           write(*,*) "initracer: error !!"
395           write(*,*) "  cannot use water option without ",
396     &                "an h2o_ice tracer !"
397           stop
398         endif
399       endif
400
401
402c     Output for records:
403c     ~~~~~~~~~~~~~~~~~~
404      write(*,*)
405      Write(*,*) '******** initracer : dust transport parameters :'
406      write(*,*) 'alpha_lift = ', alpha_lift
407      write(*,*) 'alpha_devil = ', alpha_devil
408      write(*,*) 'radius  = ', radius
409      write(*,*) 'Qext  = ', qext
410      write(*,*)
411
412      end
Note: See TracBrowser for help on using the repository browser.