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

Last change on this file since 1035 was 655, checked in by flefevre, 13 years ago
  • aeronomars/inichim_newstart.F : initialization of chemistry now handles

nitrogen species and ions (FGG, FL)

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