source: trunk/LMDZ.MARS/libf/aeronomars/inichim_newstart.F90 @ 1198

Last change on this file since 1198 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 20.3 KB
Line 
1      subroutine inichim_newstart(ngrid, nq, pq, qsurf, ps, &
2                                  flagh2o, flagthermo)
3
4      use tracer_mod
5      implicit none
6
7!=======================================================================
8!
9!  Purpose:
10!  --------
11!
12!  Initialization of the chemistry for use in the building of a new start file
13!  used by program newstart.F
14!  also used by program testphys1d.F
15!
16!  Authors:
17!  -------
18!  Initial version 11/2002 by Sebastien Lebonnois
19!  Revised 07/2003 by Monica Angelats-i-Coll to use input files
20!  Modified 10/2008 Identify tracers by their names Ehouarn Millour
21!  Modified 11/2011 Addition of methane Franck Lefevre
22!  Rewritten 04/2012 Franck Lefevre
23!
24!  Arguments:
25!  ----------
26!
27!  pq(iip1,jjp1,llm,nq)  Advected fields, ie chemical species here
28!  qsurf(ngrid,nq)     Amount of tracer on the surface (kg/m2)
29!  ps(iip1,jjp1)           Surface pressure (Pa)
30!  flagh2o                 flag for initialisation of h2o (1: yes / 0: no)
31!  flagthermo              flag for initialisation of thermosphere only (1: yes / 0: no)
32!
33!=======================================================================
34
35#include "dimensions.h"
36!#include "dimphys.h"
37#include "paramet.h"
38!#include "tracer.h"
39#include "comvert.h"
40#include "callkeys.h"
41#include "datafile.h"
42
43! inputs :
44
45      integer,intent(in) :: ngrid         ! number of atmospheric columns in the physics
46      integer,intent(in) :: nq                    ! number of tracers
47      real,intent(in) :: ps(iip1,jjp1)            ! surface pressure in the gcm (Pa)   
48      integer,intent(in) :: flagh2o               ! flag for h2o initialisation
49      integer,intent(in) :: flagthermo            ! flag for thermosphere initialisation only
50
51! outputs :
52
53      real,intent(out) :: pq(iip1,jjp1,llm,nq)  ! advected fields, ie chemical species
54      real,intent(out) :: qsurf(ngrid,nq)     ! surface values (kg/m2) of tracers
55
56! local :
57
58      integer :: iq, i, j, l, n, nbqchem
59      integer :: count, ierr, dummy
60      real    :: mmean(iip1,jjp1,llm)             ! mean molecular mass (g)
61      real    :: pgcm                             ! pressure at each layer in the gcm (Pa)
62
63      integer, parameter         :: nalt = 252    ! number of altitudes in the initialization files
64      integer                    :: nspe          ! number of species in the initialization files
65      integer, allocatable       :: niq(:)        ! local index of species in initialization files
66      real, dimension(nalt)      :: tinit, zzfile ! temperature in initialization files
67      real, dimension(nalt)      :: pinit         ! pressure in initialization files
68      real, dimension(nalt)      :: densinit      ! total number density in initialization files
69      real, allocatable          :: vmrinit(:,:)  ! mixing ratios in initialization files
70      real, allocatable          :: vmrint(:)     ! mixing ratio interpolated onto the gcm vertical grid
71      real                       :: vmr
72
73      character(len=20)          :: txt           ! to store some text
74      logical                    :: flagnitro     ! checks if N species present
75
76! 1. identify tracers by their names: (and set corresponding values of mmol)
77
78! 1.1 initialize tracer indexes to zero:
79      nqmx=nq ! initialize value of nqmx
80     
81      do iq = 1,nqmx
82        igcm_dustbin(iq) = 0
83      end do
84
85      igcm_dust_mass   = 0
86      igcm_dust_number = 0
87      igcm_ccn_mass    = 0
88      igcm_ccn_number  = 0
89      igcm_h2o_vap     = 0
90      igcm_h2o_ice     = 0
91      igcm_co2         = 0
92      igcm_co          = 0
93      igcm_o           = 0
94      igcm_o1d         = 0
95      igcm_o2          = 0
96      igcm_o3          = 0
97      igcm_h           = 0
98      igcm_h2          = 0
99      igcm_oh          = 0
100      igcm_ho2         = 0
101      igcm_h2o2        = 0
102      igcm_ch4         = 0
103      igcm_n2          = 0
104      igcm_ar          = 0
105      igcm_ar_n2       = 0
106      igcm_co2plus     = 0
107      igcm_oplus       = 0
108      igcm_o2plus      = 0
109      igcm_coplus      = 0
110      igcm_cplus       = 0
111      igcm_nplus       = 0
112      igcm_noplus      = 0
113      igcm_n2plus      = 0
114      igcm_hplus       = 0
115      igcm_hco2plus    = 0
116      igcm_elec        = 0
117
118! 1.2 find dust tracers
119
120      count = 0
121
122      if (dustbin > 0) then
123         do iq = 1,nqmx
124            txt = " "
125            write(txt,'(a4,i2.2)') 'dust', count + 1
126            if (noms(iq) == txt) then
127               count = count + 1
128               igcm_dustbin(count) = iq
129               mmol(iq) = 100.
130            end if
131         end do !do iq=1,nqmx
132      end if ! of if (dustbin.gt.0)
133
134      if (doubleq) then
135         do iq = 1,nqmx
136            if (noms(iq) == "dust_mass") then
137               igcm_dust_mass = iq
138               count = count + 1
139            end if
140            if (noms(iq) == "dust_number") then
141               igcm_dust_number = iq
142               count = count + 1
143            end if
144         end do
145      end if ! of if (doubleq)
146
147      if (scavenging) then
148         do iq = 1,nqmx
149            if (noms(iq) == "ccn_mass") then
150               igcm_ccn_mass = iq
151               count = count + 1
152            end if
153            if (noms(iq) == "ccn_number") then
154               igcm_ccn_number = iq
155               count = count + 1
156            end if
157         end do
158      end if ! of if (scavenging)
159
160      if (submicron) then
161         do iq=1,nqmx
162            if (noms(iq) == "dust_submicron") then
163               igcm_dust_submicron = iq
164               mmol(iq) = 100.
165               count = count + 1
166            end if
167         end do
168      end if ! of if (submicron)
169
170! 1.3 find chemistry and water tracers
171
172      nbqchem = 0
173      do iq = 1,nqmx
174         if (noms(iq) == "co2") then
175            igcm_co2 = iq
176            mmol(igcm_co2) = 44.
177            count = count + 1
178            nbqchem = nbqchem + 1
179        end if
180        if (noms(iq) == "co") then
181           igcm_co = iq
182           mmol(igcm_co) = 28.
183           count = count + 1
184           nbqchem = nbqchem + 1
185        end if
186        if (noms(iq) == "o") then
187           igcm_o = iq
188           mmol(igcm_o) = 16.
189           count = count + 1
190           nbqchem = nbqchem + 1
191        end if
192        if (noms(iq) == "o1d") then
193           igcm_o1d = iq
194           mmol(igcm_o1d) = 16.
195           count = count + 1
196           nbqchem = nbqchem + 1
197        end if
198        if (noms(iq) == "o2") then
199           igcm_o2 = iq
200           mmol(igcm_o2) = 32.
201           count = count + 1
202           nbqchem = nbqchem + 1
203        end if
204        if (noms(iq) == "o3") then
205           igcm_o3 = iq
206           mmol(igcm_o3) = 48.
207           count = count + 1
208           nbqchem = nbqchem + 1
209        end if
210        if (noms(iq) == "h") then
211           igcm_h = iq
212           mmol(igcm_h) = 1.
213           count = count + 1
214           nbqchem = nbqchem + 1
215        end if
216        if (noms(iq) == "h2") then
217           igcm_h2 = iq
218           mmol(igcm_h2) = 2.
219           count = count + 1
220           nbqchem = nbqchem + 1
221        end if
222        if (noms(iq) == "oh") then
223           igcm_oh = iq
224           mmol(igcm_oh) = 17.
225           count = count + 1
226           nbqchem = nbqchem + 1
227        end if
228        if (noms(iq) == "ho2") then
229           igcm_ho2 = iq
230           mmol(igcm_ho2) = 33.
231           count = count + 1
232           nbqchem = nbqchem + 1
233        end if
234        if (noms(iq) == "h2o2") then
235           igcm_h2o2 = iq
236           mmol(igcm_h2o2) = 34.
237           count = count + 1
238           nbqchem = nbqchem + 1
239        end if
240        if (noms(iq) == "ch4") then
241           igcm_ch4 = iq
242           mmol(igcm_ch4) = 16.
243           count = count + 1
244           nbqchem = nbqchem + 1
245        end if
246        if (noms(iq) == "n2") then
247           igcm_n2 = iq
248           mmol(igcm_n2) = 28.
249           count = count + 1
250           nbqchem = nbqchem + 1
251        end if
252        if (noms(iq) == "n") then
253           igcm_n = iq
254           mmol(igcm_n) = 14.
255           count = count + 1
256           nbqchem = nbqchem + 1
257        end if
258        if (noms(iq) == "n2d") then
259           igcm_n2d = iq
260           mmol(igcm_n2d) = 14.
261           count = count + 1
262           nbqchem = nbqchem + 1
263        end if
264        if (noms(iq) == "no") then
265           igcm_no = iq
266           mmol(igcm_no) = 30.
267           count = count + 1
268           nbqchem = nbqchem + 1
269        end if
270        if (noms(iq) == "no2") then
271           igcm_no2 = iq
272           mmol(igcm_no2) = 46.
273           count = count + 1
274           nbqchem = nbqchem + 1
275        end if
276        if (noms(iq) == "ar") then
277           igcm_ar = iq
278           mmol(igcm_ar) = 40.
279           count = count + 1
280           nbqchem = nbqchem + 1
281        end if
282        if (noms(iq) == "h2o_vap") then
283           igcm_h2o_vap = iq
284           mmol(igcm_h2o_vap) = 18.
285           count = count + 1
286           nbqchem = nbqchem + 1
287        end if
288        if (noms(iq) == "h2o_ice") then
289           igcm_h2o_ice = iq
290           mmol(igcm_h2o_ice) = 18.
291           count = count + 1
292           nbqchem = nbqchem + 1
293        end if
294
295! 1.4 find ions
296
297        if (noms(iq) == "co2plus") then
298           igcm_co2plus = iq
299           mmol(igcm_co2plus) = 44.
300           count = count + 1
301           nbqchem = nbqchem + 1
302        end if
303        if (noms(iq) == "oplus") then
304           igcm_oplus = iq
305           mmol(igcm_oplus) = 16.
306           count = count + 1
307           nbqchem = nbqchem + 1
308        end if
309        if (noms(iq) == "o2plus") then
310           igcm_o2plus = iq
311           mmol(igcm_o2plus) = 32.
312           count = count + 1
313           nbqchem = nbqchem + 1
314        end if
315        if (noms(iq) == "coplus") then
316           igcm_coplus = iq
317           mmol(igcm_coplus) = 28.
318           count = count + 1
319           nbqchem = nbqchem + 1
320        end if
321        if (noms(iq) == "cplus") then
322           igcm_cplus = iq
323           mmol(igcm_cplus) = 12.
324           count = count + 1
325           nbqchem = nbqchem + 1
326        end if
327        if (noms(iq) == "nplus") then
328           igcm_nplus = iq
329           mmol(igcm_nplus) = 14.
330           count = count + 1
331           nbqchem = nbqchem + 1
332        end if
333        if (noms(iq) == "noplus") then
334           igcm_noplus = iq
335           mmol(igcm_noplus) = 30.
336           count = count + 1
337           nbqchem = nbqchem + 1
338        end if
339        if (noms(iq) == "n2plus") then
340           igcm_n2plus = iq
341           mmol(igcm_n2plus) = 28.
342           count = count + 1
343           nbqchem = nbqchem + 1
344        end if
345        if (noms(iq) == "hplus") then
346           igcm_hplus = iq
347           mmol(igcm_hplus) = 1.
348           count = count + 1
349           nbqchem = nbqchem + 1
350        end if
351        if (noms(iq) == "hco2plus") then
352           igcm_hco2plus = iq
353           mmol(igcm_hco2plus) = 45.
354           count = count + 1
355           nbqchem = nbqchem + 1
356        end if
357        if (noms(iq) == "elec") then
358           igcm_elec = iq
359           mmol(igcm_elec) = 1./1822.89
360           count = count + 1
361        end if
362
363! 1.5 find idealized non-condensible tracer
364
365        if (noms(iq) == "Ar_N2") then
366           igcm_ar_n2 = iq
367           mmol(igcm_ar_n2) = 30.
368           count = count + 1
369        end if
370
371      end do ! of do iq=1,nqmx
372     
373! 1.6 check that we identified all tracers:
374
375      if (count /= nqmx) then
376         write(*,*) "inichim_newstart: found only ",count," tracers"
377         write(*,*) "                  expected ",nqmx
378         do iq = 1,count
379            write(*,*) '      ', iq, ' ', trim(noms(iq))
380         end do
381         stop
382      else
383         write(*,*) "inichim_newstart: found all expected tracers"
384         do iq = 1,nqmx
385            write(*,*) '      ', iq, ' ', trim(noms(iq))
386         end do
387      end if
388
389! 1.7 check if nitrogen species are present:
390
391      if(igcm_no == 0) then
392         !check that no N species is in traceur.def
393         if(igcm_n /= 0 .or. igcm_no2 /= 0 .or. igcm_n2d /= 0) then
394            write(*,*)'inichim_newstart error:'
395            write(*,*)'N, NO2 and/or N2D are in traceur.def, but not NO'
396            write(*,*)'stop'
397            stop
398         endif
399         flagnitro = .false.
400         nspe = 14
401      else
402         !check that all N species are in traceur.def
403         if(igcm_n == 0 .or. igcm_no2 == 0 .or. igcm_n2d == 0) then
404            write(*,*)'inichim_newstart error:'
405            write(*,*)'if NO is in traceur.def, N, NO2 and N2D must also be'
406            write(*,*)'stop'
407            stop
408         endif
409         flagnitro = .true.
410         nspe = 18
411      endif
412
413! 1.8 allocate arrays
414
415      allocate(niq(nspe))
416      allocate(vmrinit(nalt,nspe))
417      allocate(vmrint(nspe))
418
419! 2. load in chemistry data for initialization:
420
421! order of major species in initialization file:
422!
423!    1: co2
424!    2: ar
425!    3: n2 
426!    4: o2 
427!    5: co 
428!    6: o   
429!    7: h2
430!
431! order of minor species in initialization file:
432!
433!    1: h 
434!    2: oh
435!    3: ho2
436!    4: h2o
437!    5: h2o2
438!    6: o1d
439!    7: o3
440!
441! order of nitrogen species in initialization file:
442!
443!    1: n
444!    2: no
445!    3: no2
446!    4: n2d
447
448! major species:
449
450      niq(1) = igcm_co2
451      niq(2) = igcm_ar
452      niq(3) = igcm_n2
453      niq(4) = igcm_o2
454      niq(5) = igcm_co
455      niq(6) = igcm_o
456      niq(7) = igcm_h2
457
458! minor species:
459
460      niq(8)  = igcm_h
461      niq(9)  = igcm_oh
462      niq(10) = igcm_ho2
463      niq(11) = igcm_h2o_vap
464      niq(12) = igcm_h2o2
465      niq(13) = igcm_o1d
466      niq(14) = igcm_o3
467
468! nitrogen species:
469
470      if (flagnitro) then
471         niq(15) = igcm_n
472         niq(16) = igcm_no
473         niq(17) = igcm_no2
474         niq(18) = igcm_n2d         
475      end if
476
477! 2.1 open initialization files
478
479      open(210, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_may.dat')
480      if (ierr /= 0) then
481         write(*,*)'Error : cannot open file atmosfera_LMD_may.dat '
482         write(*,*)'(in aeronomars/inichim_newstart.F)'
483         write(*,*)'It should be in :', trim(datafile),'/'
484         write(*,*)'1) You can change this path in callphys.def with'
485         write(*,*)'   datadir=/path/to/datafiles/'
486         write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)'
487         write(*,*)'   can be obtained online on:'
488         write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
489         stop
490      end if
491      open(220, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_min.dat')
492      if (ierr /= 0) then
493         write(*,*)'Error : cannot open file atmosfera_LMD_min.dat '
494         write(*,*)'(in aeronomars/inichim_newstart.F)'
495         write(*,*)'It should be in :', trim(datafile),'/'
496         write(*,*)'1) You can change this path in callphys.def with'
497         write(*,*)'   datadir=/path/to/datafiles/'
498         write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)'
499         write(*,*)'   can be obtained online on:'
500         write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
501         stop
502      end if
503      if(flagnitro) then
504         open(230, iostat=ierr,file=trim(datafile)//'/atmosfera_LMD_nitr.dat')
505         if (ierr.ne.0) then
506            write(*,*)'Error : cannot open file atmosfera_LMD_nitr.dat '
507            write(*,*)'(in aeronomars/inichim_newstart.F)'
508            write(*,*)'It should be in :', datafile
509            write(*,*)'1) You can change this directory address in '
510            write(*,*)'   file phymars/datafile.h'
511            write(*,*)'2) If necessary atmosfera_LMD_nitr.dat (and others)'
512            write(*,*)'   can be obtained online on:'
513            write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
514            STOP
515         endif
516      endif   ! Of if(flagnitro)
517
518! 2.2 read initialization files
519
520! major species
521
522      read(210,*)
523      do l = 1,nalt
524         read(210,*) dummy, tinit(l), pinit(l), densinit(l), &
525                     (vmrinit(l,n), n = 1,7)
526         pinit(l) = pinit(l)*100.              ! conversion in Pa
527         pinit(l) = log(pinit(l))              ! for the vertical interpolation
528      end do
529      close(210)
530
531! minor species
532
533      read(220,*)
534      do l = 1,nalt
535         read(220,*) dummy, (vmrinit(l,n), n = 8,14)
536      end do
537      close(220)
538
539! nitrogen species
540
541      if (flagnitro) then
542         read(230,*)
543         do l = 1,nalt
544            read(230,*) dummy, (vmrinit(l,n), n = 15,18)
545         end do
546         close(230)
547      end if
548     
549! 3. initialization of tracers
550
551      do i = 1,iip1
552         do j = 1,jjp1
553            do l = 1,llm
554
555               pgcm = aps(l) + bps(l)*ps(i,j)  ! gcm pressure
556               pgcm = log(pgcm)                ! for the vertical interpolation
557               mmean(i,j,l) = 0.
558
559! 3.1 vertical interpolation
560
561               do n = 1,nspe
562                  call intrplf(pgcm,vmr,pinit,vmrinit(:,n),nalt)
563                  vmrint(n) = vmr
564                  iq = niq(n)
565                  mmean(i,j,l) = mmean(i,j,l) + vmrint(n)*mmol(iq)
566               end do
567
568! 3.2 attribute mixing ratio: - all layers or only thermosphere
569!                             - with our without h2o
570
571               if (flagthermo == 0 .or. (flagthermo == 1 .and. exp(pgcm) < 0.1)) then
572                  do n = 1,nspe
573                     iq = niq(n)
574                     if (iq /= igcm_h2o_vap .or. flagh2o == 1) then
575                        pq(i,j,l,iq) = vmrint(n)*mmol(iq)/mmean(i,j,l)
576                     end if
577                  end do
578               end if
579
580            end do
581         end do
582      end do
583
584! set surface values of chemistry tracers to zero
585
586      if (flagthermo == 0) then
587         ! NB: no problem for "surface water vapour" tracer which is always 0
588         do n = 1,nspe
589            iq = niq(n)
590            qsurf(1:ngrid,iq) = 0.
591         end do
592      end if
593
594! 3.3 initialization of tracers not contained in the initialization files
595
596! methane : 10 ppbv
597
598      if (igcm_ch4 /= 0) then
599         vmr = 10.e-9       
600         do i = 1,iip1
601            do j = 1,jjp1
602               do l = 1,llm
603                  pq(i,j,l,igcm_ch4) = vmr*mmol(igcm_ch4)/mmean(i,j,l)
604               end do
605            end do
606         end do
607         ! set surface value to zero
608         qsurf(1:ngrid,igcm_ch4) = 0.
609      end if
610
611! ions: 0
612
613      if (igcm_co2plus /= 0) then
614         !check that all required ions are in traceur.def
615         if (igcm_o2plus == 0 .or. igcm_oplus == 0 .or. igcm_coplus == 0          &
616              .or. igcm_cplus == 0 .or. igcm_nplus == 0 .or. igcm_noplus == 0    &
617              .or. igcm_n2plus == 0 .or. igcm_hplus == 0 .or. igcm_hco2plus == 0 &
618              .or. igcm_elec == 0) then
619            write(*,*)'inichim_newstart error:'
620            write(*,*)'if co2plus is in traceur.def, all other ions must also be'
621            write(*,*)'o2plus, oplus, coplus, cplus, nplus, noplus, n2plus'
622            write(*,*)'hplus, hco2plus and elec'
623            write(*,*)'stop'
624            stop
625         end if
626
627         do i = 1,iip1
628            do j = 1,jjp1
629               do l = 1,llm
630                  ! all ions to 0     
631                  pq(i,j,l,igcm_co2plus)  = 0.
632                  pq(i,j,l,igcm_o2plus)   = 0.
633                  pq(i,j,l,igcm_oplus)    = 0.
634                  pq(i,j,l,igcm_coplus)   = 0.
635                  pq(i,j,l,igcm_cplus)    = 0.
636                  pq(i,j,l,igcm_nplus)    = 0.
637                  pq(i,j,l,igcm_noplus)   = 0.
638                  pq(i,j,l,igcm_n2plus)   = 0.
639                  pq(i,j,l,igcm_hplus)    = 0.
640                  pq(i,j,l,igcm_hco2plus) = 0.
641                  pq(i,j,l,igcm_elec)     = 0.
642               end do
643            end do
644         end do
645
646         ! surface value to 0
647
648         qsurf(1:ngrid,igcm_co2plus)  = 0.
649         qsurf(1:ngrid,igcm_o2plus)   = 0.
650         qsurf(1:ngrid,igcm_oplus)    = 0.
651         qsurf(1:ngrid,igcm_coplus)   = 0.
652         qsurf(1:ngrid,igcm_cplus)    = 0.
653         qsurf(1:ngrid,igcm_nplus)    = 0.
654         qsurf(1:ngrid,igcm_noplus)   = 0.
655         qsurf(1:ngrid,igcm_n2plus)   = 0.
656         qsurf(1:ngrid,igcm_hplus)    = 0.
657         qsurf(1:ngrid,igcm_hco2plus) = 0.
658         qsurf(1:ngrid,igcm_elec)     = 0.
659
660      else
661
662         if (igcm_o2plus /= 0 .or. igcm_oplus /= 0 .or. igcm_coplus /= 0          &
663              .or. igcm_cplus /= 0 .or. igcm_nplus /= 0 .or. igcm_noplus /= 0    &
664              .or. igcm_n2plus /= 0 .or. igcm_hplus /= 0 .or. igcm_hco2plus /= 0 &
665              .or. igcm_elec /= 0) then
666            write(*,*)'inichim_newstart error:'
667            write(*,*)'some ions are in traceur.def, but not co2plus'
668            write(*,*)'stop'
669            stop
670         end if
671      end if    ! of if(igcm_co2 /= 0)
672     
673      ! deallocations
674
675      deallocate(niq)
676      deallocate(vmrinit)
677      deallocate(vmrint)
678
679      end
Note: See TracBrowser for help on using the repository browser.