source: trunk/LMDZ.PLUTO/libf/phypluto/initracer.F90 @ 3586

Last change on this file since 3586 was 3585, checked in by debatzbr, 3 weeks ago

Connecting microphysics to radiative transfer + miscellaneous cleans

File size: 20.2 KB
RevLine 
[3184]1      SUBROUTINE initracer(ngrid,nq)
2
3      use surfdat_h, ONLY: dryness
4      USE tracer_h
[3572]5      USE callkeys_mod, only: optichaze,nb_monomer,haze,fractal,fasthaze,rad_haze,callmufi
[3184]6      USE recombin_corrk_mod, ONLY: ini_recombin
7      USE mod_phys_lmdz_para, only: is_master, bcast
[3195]8      use aerosol_mod, only: iaero_haze,i_haze
[3184]9      IMPLICIT NONE
10!=======================================================================
11!   subject:
12!   --------
[3195]13!   Initialization related to tracer
[3184]14!   (transported dust, water, chemical species, ice...)
15!
16!   Name of the tracer
17!
18!   Test of dimension :
19!   Initialize COMMON tracer in tracer.h, using tracer names provided
20!   by the argument nametrac
21!
22!   author: F.Forget
23!   ------
[3557]24!            Ehouarn Millour (oct. 2008): identify tracers by their names
25!            Y. Jaziri & J. Vatant d'Ollone (2020) : modern traceur.def
26!            B. de Batz de Trenquelléon (2024): specific microphysical tracers   
[3184]27!=======================================================================
28
29      integer,intent(in) :: ngrid,nq
30
31      character(len=500) :: tracline   ! to read traceur.def lines
32      integer :: blank      !to store the index of 1st blank when reading tracers names
33      integer iq,ig,count,ierr
[3195]34      real r0_lift , reff_lift, rho_haze
35      integer nqhaze(nq)               ! to store haze tracers
[3557]36      integer i, ia, block, j
[3195]37      character(len=20) :: txt ! to store some text
[3557]38      character(LEN=20) :: tracername ! to temporarily store text
39      character(LEN=20) :: str
[3184]40
41!-----------------------------------------------------------------------
42!  radius(nq)      ! aerosol particle radius (m)
43!  rho_q(nq)       ! tracer densities (kg.m-3)
44!  qext(nq)        ! Single Scat. Extinction coeff at 0.67 um
45!  alpha_lift(nq)  ! saltation vertical flux/horiz flux ratio (m-1)
46!  alpha_devil(nq) ! lifting coeeficient by dust devil
47!  rho_dust          ! Mars dust density
48!  rho_ice           ! Water ice density
49!  doubleq           ! if method with mass (iq=1) and number(iq=2) mixing ratio
50!  varian            ! Characteristic variance of log-normal distribution
51!-----------------------------------------------------------------------
52
53      if (is_master) then ! only the master proc/thread needs do this
54
55        moderntracdef=.false. ! For modern traceur.def (default false, old type)
56
57        open(407, form = 'formatted', status = 'old', &
58             file = 'traceur.def', iostat=ierr)
59        if (ierr /=0) then
60          ! call abort_physic('initracer',&
61          ! 'Problem in opening traceur.def',1)
62          print*,'Problem in opening traceur.def'
63          stop
64        end if
65!! - Modif. by JVO and YJ for modern planetary traceur.def ---------------
66        READ(407,'(A)') tracline
67        IF (trim(tracline).ne.'#ModernTrac-v1') THEN ! Test modern traceur.def
68          READ(tracline,*) nqtot ! Try standard traceur.def
69        ELSE
70         moderntracdef = .true.
71         DO
72           READ(407,'(A)',iostat=ierr) tracline
73           IF (ierr==0) THEN
74             IF (index(tracline,'#').ne.1) THEN ! Allows arbitary number of comments lines in the header
75               READ(tracline,*) nqtot
76               ! Temporary not implemented solution
77               if (nqtot/=nq) then
78      !            call abort_physic('initracer','Different number of '// &
79      ! 'tracers in dynamics and physics not managed yet',1)
80                print*,'!= nbr oftracers in dynamics and physics not managed yet'
81                stop
82               endif
83               EXIT
84             ENDIF
85           ELSE ! If pb, or if reached EOF without having found nqtot
86      !        call abort_physic('initracer','Unable to read numbers '// &
87      !  'of tracers in traceur.def',1)
88            print*,"unable to read numbers of tracer in tracer.def"
89            stop
90           ENDIF
91         ENDDO
92        ENDIF ! if modern or standard traceur.def
[3195]93
[3184]94       endif ! of if (is_master)
[3195]95
[3184]96       ! share the information with other procs/threads (if any)
97       CALL bcast(nqtot)
98       CALL bcast(moderntracdef)
[3195]99
[3184]100!! -----------------------------------------------------------------------
101       !! For the moment number of tracers in dynamics and physics are equal
102       nqtot=nq
103       !! we allocate once for all arrays in common in tracer_h.F90
104       !! (supposedly those are not used before call to initracer)
105       IF (.NOT.ALLOCATED(noms))           ALLOCATE(noms(nq))
106       IF (.NOT.ALLOCATED(mmol))           ALLOCATE(mmol(nq))
107       IF (.NOT.ALLOCATED(aki))            ALLOCATE(aki(nqtot))
108       IF (.NOT.ALLOCATED(cpi))            ALLOCATE(cpi(nqtot))
109       IF (.NOT.ALLOCATED(radius))         ALLOCATE(radius(nq))
110       IF (.NOT.ALLOCATED(rho_q))          ALLOCATE(rho_q(nq))
111       IF (.NOT.ALLOCATED(qext))           ALLOCATE(qext(nq))
112       IF (.NOT.ALLOCATED(alpha_lift))     ALLOCATE(alpha_lift(nq))
113       IF (.NOT.ALLOCATED(alpha_devil))    ALLOCATE(alpha_devil(nq))
114       IF (.NOT.ALLOCATED(qextrhor))       ALLOCATE(qextrhor(nq))
[3195]115      !  IF (.NOT.ALLOCATED(igcm_dustbin))   ALLOCATE(igcm_dustbin(nq))
[3184]116       IF (.NOT.ALLOCATED(is_chim))        ALLOCATE(is_chim(nqtot))
117       IF (.NOT.ALLOCATED(is_rad))         ALLOCATE(is_rad(nqtot))
118       IF (.NOT.ALLOCATED(is_recomb))      ALLOCATE(is_recomb(nqtot))
119       IF (.NOT.ALLOCATED(is_recomb_qset)) THEN
120         ALLOCATE(is_recomb_qset(nqtot))
121       ENDIF
122       IF (.NOT.ALLOCATED(is_recomb_qotf)) THEN
123         ALLOCATE(is_recomb_qotf(nqtot))
124       ENDIF
125
126       !! initialization
127       alpha_lift(:)     = 0.
128       alpha_devil(:)    = 0.
129       mmol(:)           = 0.
130       aki(:)            = 0.
131       cpi(:)            = 0.
132       is_chim(:)        = 0
133       is_rad(:)         = 0
134       is_recomb(:)      = 0
135       is_recomb_qset(:) = 0
136       is_recomb_qotf(:) = 0
[3195]137
[3184]138       ! Added by JVO 2017 : these arrays are handled later
139       ! -> initialization is the least we can do, please !!!
140       radius(:)=0.
141       qext(:)=0.
142
143       rho_q(:) = 0. !need to be init here if we want to read it from modern traceur with get_tracdat
144
145
146! Initialization: Read tracers names from traceur.def
[3195]147        do iq=1,nq
[3184]148          if (is_master) read(407,'(A)') tracline
149          call bcast(tracline)
150          blank = index(tracline,' ') ! Find position of 1st blank in tracline
151          noms(iq) = tracline(1:blank) !ensure that in Modern-trac case, noms is actually the name and not all properties
152        enddo
153
154! Identify tracers by their names: (and set corresponding values of mmol)
155      ! 0. initialize tracer indexes to zero:
156      ! NB: igcm_* indexes are commons in 'tracer.h'
157      igcm_n2=0
[3195]158      igcm_ch4_gas=0
159      igcm_ch4_ice=0
160      igcm_prec_haze=0
161      igcm_co_gas=0
162      igcm_co_ice=0
[3184]163
[3195]164      nqhaze(:)=0
165      i=0
166      DO iq=1,nq
167         txt=noms(iq)
168         IF (txt(1:4).eq."haze") THEN
169            i=i+1
170            nqhaze(i)=iq
171         ENDIF
172      ENDDO
173      if ((haze.or.fasthaze).and.i==0) then
174         print*, 'Haze active but no haze tracer in traceur.def'
175         stop
176      endif
[3184]177
[3195]178      igcm_haze=0
179      igcm_haze10=0
180      igcm_haze30=0
181      igcm_haze50=0
182      igcm_haze100=0
183
184!     Eddy diffusion tracers
185      igcm_eddy1e6=0
186      igcm_eddy1e7=0
187      igcm_eddy5e7=0
188      igcm_eddy1e8=0
189      igcm_eddy5e8=0
190      write(*,*) 'initracer: noms() ', noms
191      rho_n2=1000        ! n2 ice
192      rho_ch4_ice=520.       ! rho ch4 ice
193      rho_co_ice=520.       ! rho ch4 ice
194      rho_haze=800.     ! haze
195
[3184]196      ! 1. find dust tracers
197      count=0
198
199      ! 2. find chemistry and water tracers
200      do iq=1,nq
201        if (noms(iq).eq."n2") then
202          igcm_n2=iq
203          mmol(igcm_n2)=28.
204          count=count+1
[3195]205          write(*,*) 'Tracer ',count,' = n2'
[3184]206        endif
[3275]207        if (noms(iq).eq."ch4_gas") then
[3195]208          igcm_ch4_gas=iq
209          mmol(igcm_ch4_gas)=16.
[3184]210          count=count+1
[3195]211          write(*,*) 'Tracer ',count,' = ch4 gas'
[3184]212        endif
[3195]213        if (noms(iq).eq."ch4_ice") then
214          igcm_ch4_ice=iq
215          mmol(igcm_ch4_ice)=16.
216          radius(iq)=3.e-6
217          rho_q(iq)=rho_ch4_ice
[3184]218          count=count+1
[3195]219          write(*,*) 'Tracer ',count,' = ch4 ice'
[3184]220        endif
[3275]221        if (noms(iq).eq."co_gas") then
[3195]222          igcm_co_gas=iq
223          mmol(igcm_co_gas)=28.
[3184]224          count=count+1
[3195]225          write(*,*) 'Tracer ',count,' = co gas'
[3184]226        endif
[3195]227        if (noms(iq).eq."co_ice") then
228          igcm_co_ice=iq
229          mmol(igcm_co_ice)=28.
230          radius(iq)=3.e-6
231          rho_q(iq)=rho_co_ice
[3184]232          count=count+1
[3195]233          write(*,*) 'Tracer ',count,' = co ice'
[3184]234        endif
[3195]235        if (noms(iq).eq."prec_haze") then
236          igcm_prec_haze=iq
[3184]237          count=count+1
[3195]238          write(*,*) 'Tracer ',count,' = prec haze'
[3184]239        endif
[3195]240        if (noms(iq).eq."haze") then
241          igcm_haze=iq
[3184]242          count=count+1
[3195]243          radius(iq)=rad_haze
244          rho_q(iq)=rho_haze
245          write(*,*) 'Tracer ',count,' = haze'
[3184]246        endif
[3195]247        if (noms(iq).eq."haze10") then
248          igcm_haze10=iq
[3184]249          count=count+1
[3195]250          radius(iq)=10.e-9
251          rho_q(iq)=rho_haze
252          write(*,*) 'Tracer ',count,' = haze10'
[3184]253        endif
[3195]254        if (noms(iq).eq."haze30") then
255          igcm_haze30=iq
[3184]256          count=count+1
[3195]257          radius(iq)=30.e-9
258          rho_q(iq)=rho_haze
259          write(*,*) 'Tracer ',count,' = haze30'
[3184]260        endif
[3195]261        if (noms(iq).eq."haze50") then
262          igcm_haze50=iq
[3184]263          count=count+1
[3195]264          radius(iq)=50.e-9
265          rho_q(iq)=rho_haze
266          write(*,*) 'Tracer ',count,' = haze50'
[3184]267        endif
[3195]268        if (noms(iq).eq."haze100") then
269          igcm_haze100=iq
[3184]270          count=count+1
[3195]271          radius(iq)=100.e-9
272          rho_q(iq)=rho_haze
273          write(*,*) 'Tracer ',count,' = haze100'
[3184]274        endif
[3195]275!       Eddy diffusion tracers
276        if (noms(iq).eq."eddy1e6") then
277          igcm_eddy1e6=iq
[3184]278          count=count+1
[3195]279          write(*,*) 'Tracer ',count,' = eddy1e6'
[3184]280        endif
[3195]281        if (noms(iq).eq."eddy1e7") then
282          igcm_eddy1e7=iq
[3184]283          count=count+1
[3195]284          write(*,*) 'Tracer ',count,' = eddy1e7'
[3184]285        endif
[3195]286        if (noms(iq).eq."eddy5e7") then
287          igcm_eddy5e7=iq
[3184]288          count=count+1
[3195]289          write(*,*) 'Tracer ',count,' = eddy5e7'
[3184]290        endif
[3195]291        if (noms(iq).eq."eddy1e8") then
292          igcm_eddy1e8=iq
[3184]293          count=count+1
[3195]294          write(*,*) 'Tracer ',count,' = eddy1e8'
[3184]295        endif
[3195]296        if (noms(iq).eq."eddy5e8") then
297          igcm_eddy5e8=iq
[3184]298          count=count+1
[3195]299          write(*,*) 'Tracer ',count,' = eddy5e8'
[3184]300        endif
301      enddo ! of do iq=1,nq
302
[3241]303      ! ! 3. find condensable traceurs different from h2o and n2
304      ! do iq=1,nq
305      !   if ((index(noms(iq),"vap") .ne. 0) .and. (index(noms(iq),"n2") .eq. 0)) then
306      !     count=count+1
307      !   endif
308      !   if ((index(noms(iq),"ice") .ne. 0) .and. (index(noms(iq),"n2") .eq. 0)) then
309      !     count=count+1
310      !   endif
[3184]311
[3241]312      ! enddo ! of do iq=1,nq
[3195]313
[3184]314      ! check that we identified all tracers:
315      if (count.ne.nq) then
316        write(*,*) "initracer: found only ",count," tracers"
317        write(*,*) "               expected ",nq
318        do iq=1,count
319          write(*,*)'      ',iq,' ',trim(noms(iq))
320        enddo
321!        stop
322      else
323        write(*,*) "initracer: found all expected tracers, namely:"
324        do iq=1,nq
325          write(*,*)'      ',iq,' ',trim(noms(iq))
326        enddo
327      endif
[3557]328     
329      ! Compute number of microphysics tracers:
330      ! By convention they all have the prefix "mu_" (case sensitive !)
331      nmicro = 0
332      IF (callmufi) THEN
333         DO iq=1,nq
334            str = noms(iq)
335            IF (str(1:3) == "mu_") THEN
336               nmicro = nmicro+1
337               count = count+1
338            ENDIF
339         ENDDO
[3184]340
[3557]341         ! Checking the expected number of tracers:
342         ! Microphysics moment model: nmicro = 4
343         IF (nmicro < 4) THEN
344            WRITE(*,*) "initracer:error:"," Inconsistent number of microphysical tracers"
345            WRITE(*,*) "expected at least 4 tracers,", nmicro, " given"
346            CALL abort
347         ELSE IF (nmicro > 4) THEN
348            WRITE(*,*) "!!! WARNING !!! initracer: I was expecting only four tracers, you gave me more."
349            CALL abort
350         ENDIF
351
352         ! microphysics indexes share the same values than original tracname.
353         IF (.NOT.ALLOCATED(micro_indx)) ALLOCATE(micro_indx(nmicro))
354         j = 1
355         DO i=1,nq
356            str = noms(i)
357            IF (str(1:3) == "mu_") THEN
358               micro_indx(j) = i
359               j=j+1
360            ENDIF
361         ENDDO
362     
363      ELSE
364         IF (.NOT.ALLOCATED(micro_indx)) ALLOCATE(micro_indx(nmicro))
365     
366      ENDIF ! end of callmufi
367
[3184]368      ! Get data of tracers. Need to rewind traceur.def first
369      if (is_master) then
370       rewind(407)
[3195]371       do
[3184]372        read(407,'(A)') tracline
[3195]373        if (index(tracline,"#") .ne. 1) then
[3184]374          exit
375        endif
376       enddo
377      endif
378      do iq=1,nqtot
379        if (is_master) read(407,'(A)') tracline
380        call bcast(tracline)
381        call get_tracdat(iq, tracline)
382      enddo
383
384      if (is_master) close(407)
385
386      ! Calculate number of species in the chemistry
387      nesp = sum(is_chim)
388      write(*,*) 'Number of species in the chemistry nesp = ',nesp
389
[3557]390      ! Calculate number of microphysical tracer
391      write(*,*) 'Number of microphysical tracer nmicro = ',nmicro
392      IF (callmufi) THEN
393         call dumptracers(micro_indx)
394      ENDIF
395
[3184]396!     Processing modern traceur options
397      if(moderntracdef) then
398        call ini_recombin
399      endif
400
401!------------------------------------------------------------
402!     Initialisation tracers ....
403!------------------------------------------------------------
404      ! rho_q(1:nq)=0
405
406      rho_n2=1000.   ! N2 ice density (kg.m-3)
407
[3193]408      lw_co=274000.
409      lw_ch4=586700.
410      lw_n2=2.5e5
[3237]411      write(*,*) "lw_n2 = ", lw_n2
[3195]412
[3585]413      if (callmufi) then
414        if (optichaze) then
415          iaero_haze = 2
416          write(*,*) 'Microphysical moment model'
417          write(*,*) '--- number of haze aerosol = ', iaero_haze
418        endif ! end optichaze
419     
420      elseif (haze) then
[3193]421        ! the sedimentation radius remains radius(igcm_haze)
422        if (fractal) then
[3585]423          nmono=nb_monomer
[3193]424        else
[3585]425          nmono=1
[3557]426        endif ! end fractal
[3585]427       
428        ia = 0
[3572]429        if (optichaze) then
[3585]430          ia = ia + 1
431          iaero_haze = ia
[3572]432          write(*,*) '--- number of haze aerosol = ', iaero_haze
[3585]433          block = 0 ! Only one type of haze is active : the first one set in traceur.def
[3572]434          do iq=1,nq
435            tracername=noms(iq)
436            write(*,*) "--> tracername ",iq,'/',nq,' = ',tracername
437            if (tracername(1:4).eq."haze".and.block.eq.0) then
438              i_haze=iq
439              block=1
440              write(*,*) "i_haze=",i_haze
441              write(*,*) "Careful: if you set many haze traceurs in &
[3585]442              traceur.def,only ",tracername," will be radiatively active &
443              (first one in traceur.def)"
[3572]444            endif
[3585]445            enddo
[3572]446        endif ! end optichaze
447      endif ! end callmufi or haze
[3193]448
[3184]449!     Output for records:
450!     ~~~~~~~~~~~~~~~~~~
451      write(*,*)
452      Write(*,*) '******** initracer : dust transport parameters :'
453      write(*,*) 'alpha_lift = ', alpha_lift
454      write(*,*) 'alpha_devil = ', alpha_devil
455      write(*,*) 'radius  = ', radius
[3195]456      write(*,*) 'Qext  = ', qext
[3184]457      write(*,*)
458
459      contains
460
461      subroutine get_tracdat(iq,tracline)
462        !-------------------ADDING NEW OPTIONS-------------------
463        ! Duplicate if sentence for adding new options
464        ! Do not forget to add the new variables in tracer_h.F90
465        ! Do not forget to allocate and initialize the new variables above
466        ! Please update list of options in "LMDZ.GENERIC/deftank/traceur.def.modern"
467        !-------------------------------------------------------
468          use tracer_h
469          implicit none
470          integer,           intent(in) :: iq       ! tracer index
471          character(len=500),intent(in) :: tracline ! traceur.def lines with parameters
[3195]472
[3184]473          read(tracline,*) noms(iq)
474          ! JVO 20 : We should add a sanity check aborting when duplicates in names !
475          write(*,*)"initracer: iq=",iq,"noms(iq)=",trim(noms(iq))
476          ! option mmol
477          if (index(tracline,'mmol='   ) /= 0) then
478              read(tracline(index(tracline,'mmol=')+len('mmol='):),*)&
479                  mmol(iq)
480              write(*,*) ' Parameter value (traceur.def) : mmol=', &
481                  mmol(iq)
482          else
483              write(*,*) ' Parameter value (default)     : mmol=', &
484                 mmol(iq)
485          end if
[3195]486          ! option aki
[3184]487          if (index(tracline,'aki='   ) /= 0) then
488              read(tracline(index(tracline,'aki=')+len('aki='):),*) &
489                  aki(iq)
490              write(*,*) ' Parameter value (traceur.def) : aki=', &
491                  aki(iq)
492          else
493              write(*,*) ' Parameter value (default)     : aki=', &
494                  aki(iq)
495          end if
[3195]496          ! option cpi
[3184]497          if (index(tracline,'cpi='   ) /= 0) then
498              read(tracline(index(tracline,'cpi=')+len('cpi='):),*) &
499                  cpi(iq)
500              write(*,*) ' Parameter value (traceur.def) : cpi=', &
501                  cpi(iq)
502          else
503              write(*,*) ' Parameter value (default)     : cpi=', &
504                  cpi(iq)
505          end if
[3195]506          ! option is_chim
[3184]507          if (index(tracline,'is_chim='   ) /= 0) then
508          read(tracline(index(tracline,'is_chim=')+len('is_chim='):),*) &
509                  is_chim(iq)
510              write(*,*) ' Parameter value (traceur.def) : is_chim=', &
511                  is_chim(iq)
512          else
513              write(*,*) ' Parameter value (default)     : is_chim=', &
514                  is_chim(iq)
515          end if
[3195]516          ! option is_rad
[3184]517          if (index(tracline,'is_rad=') /= 0) then
518              read(tracline(index(tracline,'is_rad=') &
519              +len('is_rad='):),*) is_rad(iq)
520              write(*,*) ' Parameter value (traceur.def) : is_rad=', &
521              is_rad(iq)
522          else
523              write(*,*) ' Parameter value (default)     : is_rad=',  &
524              is_rad(iq)
525          end if
526          ! option is_recomb
527          if (index(tracline,'is_recomb=') /= 0) then
528              read(tracline(index(tracline,'is_recomb=') &
529              +len('is_recomb='):),*) is_recomb(iq)
530              write(*,*) ' Parameter value (traceur.def) : is_recomb=', &
531              is_recomb(iq)
532          else
533              write(*,*) ' Parameter value (default)     : is_recomb=', &
534              is_recomb(iq)
535          end if
536          ! option is_recomb_qset
537          if (index(tracline,'is_recomb_qset=') /= 0) then
538              read(tracline(index(tracline,'is_recomb_qset=') &
539              +len('is_recomb_qset='):),*) is_recomb_qset(iq)
540              write(*,*) ' Parameter value (traceur.def) :'// &
541              ' is_recomb_qset=', &
542              is_recomb_qset(iq)
543          else
544              write(*,*) ' Parameter value (default)     :'// &
545              ' is_recomb_qset=', &
546              is_recomb_qset(iq)
547          endif
548          ! option is_recomb_qotf
549          if (index(tracline,'is_recomb_qotf=') /= 0) then
550              read(tracline(index(tracline,'is_recomb_qotf=') &
551              +len('is_recomb_qotf='):),*) is_recomb_qotf(iq)
552              write(*,*) ' Parameter value (traceur.def) :'// &
553              ' is_recomb_qotf=', &
554              is_recomb_qotf(iq)
555          else
556              write(*,*) ' Parameter value (default)     :'// &
557              ' is_recomb_qotf=',  &
558              is_recomb_qotf(iq)
559          end if
560          !option radius
[3195]561          if (index(tracline,'radius=') .ne. 0) then
[3184]562            read(tracline(index(tracline,'radius=') &
563              +len('radius='):),*) radius(iq)
564            write(*,*)'Parameter value (traceur.def) :'// &
565              "radius=",radius(iq), "m"
566          else
567            write(*,*) ' Parameter value (default)     :'// &
[3195]568            ' radius=', radius(iq)," m"
[3184]569          endif
570          !option rho
[3195]571          if (index(tracline,'rho=') .ne. 0) then
[3184]572            read(tracline(index(tracline,'rho=') &
573              +len('rho='):),*) rho_q(iq)
574            write(*,*)'Parameter value (traceur.def) :'// &
575              "rho=",rho_q(iq)
576          else
577            write(*,*) ' Parameter value (default)     :'// &
[3195]578              ' rho=', rho_q(iq)
[3184]579          endif
580      end subroutine get_tracdat
581
582      end
583
Note: See TracBrowser for help on using the repository browser.