source: branches/FLEXPART_9.1.3/src/concoutput.f90

Last change on this file was 13, checked in by saeck, 11 years ago

update to wetdepo.f90

File size: 21.8 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine concoutput(itime,outnum,gridtotalunc,wetgridtotalunc, &
23       drygridtotalunc)
24  !                        i     i          o             o
25  !       o
26  !*****************************************************************************
27  !                                                                            *
28  !     Output of the concentration grid and the receptor concentrations.      *
29  !                                                                            *
30  !     Author: A. Stohl                                                       *
31  !                                                                            *
32  !     24 May 1995                                                            *
33  !                                                                            *
34  !     13 April 1999, Major update: if output size is smaller, dump output    *
35  !                    in sparse matrix format; additional output of           *
36  !                    uncertainty                                             *
37  !                                                                            *
38  !     05 April 2000, Major update: output of age classes; output for backward*
39  !                    runs is time spent in grid cell times total mass of     *
40  !                    species.                                                *
41  !                                                                            *
42  !     17 February 2002, Appropriate dimensions for backward and forward runs *
43  !                       are now specified in file par_mod                    *
44  !                                                                            *
45  !     June 2006, write grid in sparse matrix with a single write command     *
46  !                in order to save disk space                                 *
47  !                                                                            *
48  !     2008 new sparse matrix format                                          *
49  !                                                                            *
50  !*****************************************************************************
51  !                                                                            *
52  ! Variables:                                                                 *
53  ! outnum          number of samples                                          *
54  ! ncells          number of cells with non-zero concentrations               *
55  ! sparse          .true. if in sparse matrix format, else .false.            *
56  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
57  !                                                                            *
58  !*****************************************************************************
59
60  use unc_mod
61  use point_mod
62  use outg_mod
63  use par_mod
64  use com_mod
65
66  implicit none
67
68  real(kind=dp) :: jul
69  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
70  integer :: sp_count_i,sp_count_r
71  real :: sp_fact
72  real :: outnum,densityoutrecept(maxreceptor),xl,yl
73
74  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
75  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
76  !    +    maxageclass)
77  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
78  !    +       maxageclass)
79  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
80  !    +       maxpointspec_act,maxageclass)
81  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
82  !    +       maxpointspec_act,maxageclass),
83  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
84  !    +     maxpointspec_act,maxageclass),
85  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
86  !    +     maxpointspec_act,maxageclass)
87  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
88  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
89  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
90
91  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
92  real :: auxgrid(nclassunc),gridtotal,gridsigmatotal,gridtotalunc
93  real :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
94  real :: drygridtotal,drygridsigmatotal,drygridtotalunc
95  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
96  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
97  real,parameter :: weightair=28.97
98  logical :: sp_zer
99  character :: adate*8,atime*6
100  character(len=3) :: anspec
101
102
103  ! Determine current calendar date, needed for the file name
104  !**********************************************************
105
106  jul=bdate+real(itime,kind=dp)/86400._dp
107  call caldate(jul,jjjjmmdd,ihmmss)
108  write(adate,'(i8.8)') jjjjmmdd
109  write(atime,'(i6.6)') ihmmss
110  write(unitdates,'(a)') adate//atime
111
112
113  ! For forward simulations, output fields have dimension MAXSPEC,
114  ! for backward simulations, output fields have dimension MAXPOINT.
115  ! Thus, make loops either about nspec, or about numpoint
116  !*****************************************************************
117
118
119  if (ldirect.eq.1) then
120    do ks=1,nspec
121      do kp=1,maxpointspec_act
122        tot_mu(ks,kp)=1
123      end do
124    end do
125  else
126    do ks=1,nspec
127      do kp=1,maxpointspec_act
128        tot_mu(ks,kp)=xmass(kp,ks)
129      end do
130    end do
131  endif
132
133
134  !*******************************************************************
135  ! Compute air density: sufficiently accurate to take it
136  ! from coarse grid at some time
137  ! Determine center altitude of output layer, and interpolate density
138  ! data to that altitude
139  !*******************************************************************
140
141  do kz=1,numzgrid
142    if (kz.eq.1) then
143      halfheight=outheight(1)/2.
144    else
145      halfheight=(outheight(kz)+outheight(kz-1))/2.
146    endif
147    do kzz=2,nz
148      if ((height(kzz-1).lt.halfheight).and. &
149           (height(kzz).gt.halfheight)) goto 46
150    end do
15146   kzz=max(min(kzz,nz),2)
152    dz1=halfheight-height(kzz-1)
153    dz2=height(kzz)-halfheight
154    dz=dz1+dz2
155    do jy=0,numygrid-1
156      do ix=0,numxgrid-1
157        xl=outlon0+real(ix)*dxout
158        yl=outlat0+real(jy)*dyout
159        xl=(xl-xlon0)/dx
160        yl=(yl-ylat0)/dy !v9.1.1
161        iix=max(min(nint(xl),nxmin1),0)
162        jjy=max(min(nint(yl),nymin1),0)
163        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
164             rho(iix,jjy,kzz-1,2)*dz2)/dz
165      end do
166    end do
167  end do
168
169    do i=1,numreceptor
170      xl=xreceptor(i)
171      yl=yreceptor(i)
172      iix=max(min(nint(xl),nxmin1),0)
173      jjy=max(min(nint(yl),nymin1),0)
174      densityoutrecept(i)=rho(iix,jjy,1,2)
175    end do
176
177
178  ! Output is different for forward and backward simulations
179    do kz=1,numzgrid
180      do jy=0,numygrid-1
181        do ix=0,numxgrid-1
182          if (ldirect.eq.1) then
183            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
184          else
185            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
186          endif
187        end do
188      end do
189    end do
190
191  !*********************************************************************
192  ! Determine the standard deviation of the mean concentration or mixing
193  ! ratio (uncertainty of the output) and the dry and wet deposition
194  !*********************************************************************
195
196  gridtotal=0.
197  gridsigmatotal=0.
198  gridtotalunc=0.
199  wetgridtotal=0.
200  wetgridsigmatotal=0.
201  wetgridtotalunc=0.
202  drygridtotal=0.
203  drygridsigmatotal=0.
204  drygridtotalunc=0.
205
206  do ks=1,nspec
207
208  write(anspec,'(i3.3)') ks
209  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
210    if (ldirect.eq.1) then
211      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
212           atime//'_'//anspec,form='unformatted')
213    else
214      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
215           atime//'_'//anspec,form='unformatted')
216    endif
217     write(unitoutgrid) itime
218   endif
219
220  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
221   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
222        atime//'_'//anspec,form='unformatted')
223
224    write(unitoutgridppt) itime
225  endif
226
227  do kp=1,maxpointspec_act
228  do nage=1,nageclass
229
230    do jy=0,numygrid-1
231      do ix=0,numxgrid-1
232
233  ! WET DEPOSITION
234        if ((WETDEP).and.(ldirect.gt.0)) then
235            do l=1,nclassunc
236              auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
237            end do
238            call mean(auxgrid,wetgrid(ix,jy), &
239                 wetgridsigma(ix,jy),nclassunc)
240  ! Multiply by number of classes to get total concentration
241            wetgrid(ix,jy)=wetgrid(ix,jy) &
242                 *nclassunc
243            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
244  ! Calculate standard deviation of the mean
245            wetgridsigma(ix,jy)= &
246                 wetgridsigma(ix,jy)* &
247                 sqrt(real(nclassunc))
248            wetgridsigmatotal=wetgridsigmatotal+ &
249                 wetgridsigma(ix,jy)
250        endif
251
252  ! DRY DEPOSITION
253        if ((DRYDEP).and.(ldirect.gt.0)) then
254            do l=1,nclassunc
255              auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
256            end do
257            call mean(auxgrid,drygrid(ix,jy), &
258                 drygridsigma(ix,jy),nclassunc)
259  ! Multiply by number of classes to get total concentration
260            drygrid(ix,jy)=drygrid(ix,jy)* &
261                 nclassunc
262            drygridtotal=drygridtotal+drygrid(ix,jy)
263  ! Calculate standard deviation of the mean
264            drygridsigma(ix,jy)= &
265                 drygridsigma(ix,jy)* &
266                 sqrt(real(nclassunc))
267125         drygridsigmatotal=drygridsigmatotal+ &
268                 drygridsigma(ix,jy)
269        endif
270
271  ! CONCENTRATION OR MIXING RATIO
272        do kz=1,numzgrid
273            do l=1,nclassunc
274              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
275            end do
276            call mean(auxgrid,grid(ix,jy,kz), &
277                 gridsigma(ix,jy,kz),nclassunc)
278  ! Multiply by number of classes to get total concentration
279            grid(ix,jy,kz)= &
280                 grid(ix,jy,kz)*nclassunc
281            gridtotal=gridtotal+grid(ix,jy,kz)
282  ! Calculate standard deviation of the mean
283            gridsigma(ix,jy,kz)= &
284                 gridsigma(ix,jy,kz)* &
285                 sqrt(real(nclassunc))
286            gridsigmatotal=gridsigmatotal+ &
287                 gridsigma(ix,jy,kz)
288        end do
289      end do
290    end do
291
292
293
294
295  !*******************************************************************
296  ! Generate output: may be in concentration (ng/m3) or in mixing
297  ! ratio (ppt) or both
298  ! Output the position and the values alternated multiplied by
299  ! 1 or -1, first line is number of values, number of positions
300  ! For backward simulations, the unit is seconds, stored in grid_time
301  !*******************************************************************
302
303  ! Concentration output
304  !*********************
305  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
306
307  ! Wet deposition
308         sp_count_i=0
309         sp_count_r=0
310         sp_fact=-1.
311         sp_zer=.true.
312         if ((ldirect.eq.1).and.(WETDEP)) then
313         do jy=0,numygrid-1
314            do ix=0,numxgrid-1
315  !oncentraion greater zero
316              if (wetgrid(ix,jy).gt.smallnum) then
317                 if (sp_zer.eqv..true.) then ! first non zero value
318                    sp_count_i=sp_count_i+1
319                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
320                    sp_zer=.false.
321                    sp_fact=sp_fact*(-1.)
322                 endif
323                 sp_count_r=sp_count_r+1
324                 sparse_dump_r(sp_count_r)= &
325                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
326  !                sparse_dump_u(sp_count_r)=
327  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
328              else ! concentration is zero
329                  sp_zer=.true.
330              endif
331            end do
332         end do
333         else
334            sp_count_i=0
335            sp_count_r=0
336         endif
337         write(unitoutgrid) sp_count_i
338         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
339         write(unitoutgrid) sp_count_r
340         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
341  !       write(unitoutgrid) sp_count_u
342  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
343
344  ! Dry deposition
345         sp_count_i=0
346         sp_count_r=0
347         sp_fact=-1.
348         sp_zer=.true.
349         if ((ldirect.eq.1).and.(DRYDEP)) then
350          do jy=0,numygrid-1
351            do ix=0,numxgrid-1
352              if (drygrid(ix,jy).gt.smallnum) then
353                 if (sp_zer.eqv..true.) then ! first non zero value
354                    sp_count_i=sp_count_i+1
355                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
356                    sp_zer=.false.
357                    sp_fact=sp_fact*(-1.)
358                 endif
359                 sp_count_r=sp_count_r+1
360                 sparse_dump_r(sp_count_r)= &
361                      sp_fact* &
362                      1.e12*drygrid(ix,jy)/area(ix,jy)
363  !                sparse_dump_u(sp_count_r)=
364  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
365              else ! concentration is zero
366                  sp_zer=.true.
367              endif
368            end do
369          end do
370         else
371            sp_count_i=0
372            sp_count_r=0
373         endif
374         write(unitoutgrid) sp_count_i
375         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
376         write(unitoutgrid) sp_count_r
377         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
378  !       write(*,*) sp_count_u
379  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
380
381
382
383  ! Concentrations
384         sp_count_i=0
385         sp_count_r=0
386         sp_fact=-1.
387         sp_zer=.true.
388          do kz=1,numzgrid
389            do jy=0,numygrid-1
390              do ix=0,numxgrid-1
391                if (grid(ix,jy,kz).gt.smallnum) then
392                  if (sp_zer.eqv..true.) then ! first non zero value
393                    sp_count_i=sp_count_i+1
394                    sparse_dump_i(sp_count_i)= &
395                         ix+jy*numxgrid+kz*numxgrid*numygrid
396                    sp_zer=.false.
397                    sp_fact=sp_fact*(-1.)
398                   endif
399                   sp_count_r=sp_count_r+1
400                   sparse_dump_r(sp_count_r)= &
401                        sp_fact* &
402                        grid(ix,jy,kz)* &
403                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
404  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
405  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
406  !                sparse_dump_u(sp_count_r)=
407  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
408  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
409              else ! concentration is zero
410                  sp_zer=.true.
411              endif
412              end do
413            end do
414          end do
415         write(unitoutgrid) sp_count_i
416         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
417         write(unitoutgrid) sp_count_r
418         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
419  !       write(unitoutgrid) sp_count_u
420  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
421
422
423
424    endif !  concentration output
425
426  ! Mixing ratio output
427  !********************
428
429  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
430
431  ! Wet deposition
432         sp_count_i=0
433         sp_count_r=0
434         sp_fact=-1.
435         sp_zer=.true.
436         if ((ldirect.eq.1).and.(WETDEP)) then
437          do jy=0,numygrid-1
438            do ix=0,numxgrid-1
439                if (wetgrid(ix,jy).gt.smallnum) then
440                  if (sp_zer.eqv..true.) then ! first non zero value
441                    sp_count_i=sp_count_i+1
442                    sparse_dump_i(sp_count_i)= &
443                         ix+jy*numxgrid
444                    sp_zer=.false.
445                    sp_fact=sp_fact*(-1.)
446                 endif
447                 sp_count_r=sp_count_r+1
448                 sparse_dump_r(sp_count_r)= &
449                      sp_fact* &
450                      1.e12*wetgrid(ix,jy)/area(ix,jy)
451  !                sparse_dump_u(sp_count_r)=
452  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
453              else ! concentration is zero
454                  sp_zer=.true.
455              endif
456            end do
457          end do
458         else
459           sp_count_i=0
460           sp_count_r=0
461         endif
462         write(unitoutgridppt) sp_count_i
463         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
464         write(unitoutgridppt) sp_count_r
465         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
466  !       write(unitoutgridppt) sp_count_u
467  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
468
469
470  ! Dry deposition
471         sp_count_i=0
472         sp_count_r=0
473         sp_fact=-1.
474         sp_zer=.true.
475         if ((ldirect.eq.1).and.(DRYDEP)) then
476          do jy=0,numygrid-1
477            do ix=0,numxgrid-1
478                if (drygrid(ix,jy).gt.smallnum) then
479                  if (sp_zer.eqv..true.) then ! first non zero value
480                    sp_count_i=sp_count_i+1
481                    sparse_dump_i(sp_count_i)= &
482                         ix+jy*numxgrid
483                    sp_zer=.false.
484                    sp_fact=sp_fact*(-1)
485                 endif
486                 sp_count_r=sp_count_r+1
487                 sparse_dump_r(sp_count_r)= &
488                      sp_fact* &
489                      1.e12*drygrid(ix,jy)/area(ix,jy)
490  !                sparse_dump_u(sp_count_r)=
491  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
492              else ! concentration is zero
493                  sp_zer=.true.
494              endif
495            end do
496          end do
497         else
498           sp_count_i=0
499           sp_count_r=0
500         endif
501         write(unitoutgridppt) sp_count_i
502         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
503         write(unitoutgridppt) sp_count_r
504         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
505  !       write(unitoutgridppt) sp_count_u
506  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
507
508
509  ! Mixing ratios
510         sp_count_i=0
511         sp_count_r=0
512         sp_fact=-1.
513         sp_zer=.true.
514          do kz=1,numzgrid
515            do jy=0,numygrid-1
516              do ix=0,numxgrid-1
517                if (grid(ix,jy,kz).gt.smallnum) then
518                  if (sp_zer.eqv..true.) then ! first non zero value
519                    sp_count_i=sp_count_i+1
520                    sparse_dump_i(sp_count_i)= &
521                         ix+jy*numxgrid+kz*numxgrid*numygrid
522                    sp_zer=.false.
523                    sp_fact=sp_fact*(-1.)
524                 endif
525                 sp_count_r=sp_count_r+1
526                 sparse_dump_r(sp_count_r)= &
527                      sp_fact* &
528                      1.e12*grid(ix,jy,kz) &
529                      /volume(ix,jy,kz)/outnum* &
530                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
531  !                sparse_dump_u(sp_count_r)=
532  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
533  !+              outnum*weightair/weightmolar(ks)/
534  !+              densityoutgrid(ix,jy,kz)
535              else ! concentration is zero
536                  sp_zer=.true.
537              endif
538              end do
539            end do
540          end do
541         write(unitoutgridppt) sp_count_i
542         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
543         write(unitoutgridppt) sp_count_r
544         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
545  !       write(unitoutgridppt) sp_count_u
546  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
547
548      endif ! output for ppt
549
550  end do
551  end do
552
553    close(unitoutgridppt)
554    close(unitoutgrid)
555
556  end do
557
558  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
559  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
560       wetgridtotal
561  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
562       drygridtotal
563
564  ! Dump of receptor concentrations
565
566    if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
567      write(unitoutreceptppt) itime
568      do ks=1,nspec
569        write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
570             weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
571      end do
572    endif
573
574  ! Dump of receptor concentrations
575
576    if (numreceptor.gt.0) then
577      write(unitoutrecept) itime
578      do ks=1,nspec
579        write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
580             i=1,numreceptor)
581      end do
582    endif
583
584
585
586  ! Reinitialization of grid
587  !*************************
588
589  do ks=1,nspec
590  do kp=1,maxpointspec_act
591    do i=1,numreceptor
592      creceptor(i,ks)=0.
593    end do
594    do jy=0,numygrid-1
595      do ix=0,numxgrid-1
596        do l=1,nclassunc
597          do nage=1,nageclass
598            do kz=1,numzgrid
599              gridunc(ix,jy,kz,ks,kp,l,nage)=0.
600            end do
601          end do
602        end do
603      end do
604    end do
605  end do
606  end do
607
608
609end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG