source: branches/flexpart91_hasod/src_parallel/concoutput.f90 @ 8

Last change on this file since 8 was 8, checked in by hasod, 11 years ago

Added parallel version of Flexpart91

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  ! Determine current calendar date, needed for the file name
103  !**********************************************************
104
105  jul=bdate+real(itime,kind=dp)/86400._dp
106  call caldate(jul,jjjjmmdd,ihmmss)
107  write(adate,'(i8.8)') jjjjmmdd
108  write(atime,'(i6.6)') ihmmss
109  write(unitdates,'(a)') adate//atime
110
111
112  ! For forward simulations, output fields have dimension MAXSPEC,
113  ! for backward simulations, output fields have dimension MAXPOINT.
114  ! Thus, make loops either about nspec, or about numpoint
115  !*****************************************************************
116
117
118  if (ldirect.eq.1) then
119    do ks=1,nspec
120      do kp=1,maxpointspec_act
121        tot_mu(ks,kp)=1
122      end do
123    end do
124  else
125    do ks=1,nspec
126      do kp=1,maxpointspec_act
127        tot_mu(ks,kp)=xmass(kp,ks)
128      end do
129    end do
130  endif
131
132
133  !*******************************************************************
134  ! Compute air density: sufficiently accurate to take it
135  ! from coarse grid at some time
136  ! Determine center altitude of output layer, and interpolate density
137  ! data to that altitude
138  !*******************************************************************
139
140  do kz=1,numzgrid
141    if (kz.eq.1) then
142      halfheight=outheight(1)/2.
143    else
144      halfheight=(outheight(kz)+outheight(kz-1))/2.
145    endif
146    do kzz=2,nz
147      if ((height(kzz-1).lt.halfheight).and. &
148           (height(kzz).gt.halfheight)) goto 46
149    end do
15046   kzz=max(min(kzz,nz),2)
151    dz1=halfheight-height(kzz-1)
152    dz2=height(kzz)-halfheight
153    dz=dz1+dz2
154    do jy=0,numygrid-1
155      do ix=0,numxgrid-1
156        xl=outlon0+real(ix)*dxout
157        yl=outlat0+real(jy)*dyout
158        xl=(xl-xlon0)/dx
159        yl=(yl-ylat0)/dx
160        iix=max(min(nint(xl),nxmin1),0)
161        jjy=max(min(nint(yl),nymin1),0)
162        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
163             rho(iix,jjy,kzz-1,2)*dz2)/dz
164      end do
165    end do
166  end do
167
168    do i=1,numreceptor
169      xl=xreceptor(i)
170      yl=yreceptor(i)
171      iix=max(min(nint(xl),nxmin1),0)
172      jjy=max(min(nint(yl),nymin1),0)
173      densityoutrecept(i)=rho(iix,jjy,1,2)
174    end do
175
176
177  ! Output is different for forward and backward simulations
178    do kz=1,numzgrid
179      do jy=0,numygrid-1
180        do ix=0,numxgrid-1
181          if (ldirect.eq.1) then
182            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
183          else
184            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
185          endif
186        end do
187      end do
188    end do
189
190  !*********************************************************************
191  ! Determine the standard deviation of the mean concentration or mixing
192  ! ratio (uncertainty of the output) and the dry and wet deposition
193  !*********************************************************************
194
195  gridtotal=0.
196  gridsigmatotal=0.
197  gridtotalunc=0.
198  wetgridtotal=0.
199  wetgridsigmatotal=0.
200  wetgridtotalunc=0.
201  drygridtotal=0.
202  drygridsigmatotal=0.
203  drygridtotalunc=0.
204
205  do ks=1,nspec
206
207  write(anspec,'(i3.3)') ks
208  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
209    if (ldirect.eq.1) then
210      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
211           atime//'_'//anspec,form='unformatted')
212    else
213      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
214           atime//'_'//anspec,form='unformatted')
215    endif
216     write(unitoutgrid) itime
217   endif
218
219  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
220   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
221        atime//'_'//anspec,form='unformatted')
222
223    write(unitoutgridppt) itime
224  endif
225
226  do kp=1,maxpointspec_act
227  do nage=1,nageclass
228
229    do jy=0,numygrid-1
230      do ix=0,numxgrid-1
231
232  ! WET DEPOSITION
233        if ((WETDEP).and.(ldirect.gt.0)) then
234            do l=1,nclassunc
235              auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
236            end do
237            call mean(auxgrid,wetgrid(ix,jy), &
238                 wetgridsigma(ix,jy),nclassunc)
239  ! Multiply by number of classes to get total concentration
240            wetgrid(ix,jy)=wetgrid(ix,jy) &
241                 *nclassunc
242            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
243  ! Calculate standard deviation of the mean
244            wetgridsigma(ix,jy)= &
245                 wetgridsigma(ix,jy)* &
246                 sqrt(real(nclassunc))
247            wetgridsigmatotal=wetgridsigmatotal+ &
248                 wetgridsigma(ix,jy)
249        endif
250
251  ! DRY DEPOSITION
252        if ((DRYDEP).and.(ldirect.gt.0)) then
253            do l=1,nclassunc
254              auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
255            end do
256            call mean(auxgrid,drygrid(ix,jy), &
257                 drygridsigma(ix,jy),nclassunc)
258  ! Multiply by number of classes to get total concentration
259            drygrid(ix,jy)=drygrid(ix,jy)* &
260                 nclassunc
261            drygridtotal=drygridtotal+drygrid(ix,jy)
262  ! Calculate standard deviation of the mean
263            drygridsigma(ix,jy)= &
264                 drygridsigma(ix,jy)* &
265                 sqrt(real(nclassunc))
266125         drygridsigmatotal=drygridsigmatotal+ &
267                 drygridsigma(ix,jy)
268        endif
269
270  ! CONCENTRATION OR MIXING RATIO
271        do kz=1,numzgrid
272            do l=1,nclassunc
273              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
274            end do
275            call mean(auxgrid,grid(ix,jy,kz), &
276                 gridsigma(ix,jy,kz),nclassunc)
277  ! Multiply by number of classes to get total concentration
278            grid(ix,jy,kz)= &
279                 grid(ix,jy,kz)*nclassunc
280            gridtotal=gridtotal+grid(ix,jy,kz)
281  ! Calculate standard deviation of the mean
282            gridsigma(ix,jy,kz)= &
283                 gridsigma(ix,jy,kz)* &
284                 sqrt(real(nclassunc))
285            gridsigmatotal=gridsigmatotal+ &
286                 gridsigma(ix,jy,kz)
287        end do
288      end do
289    end do
290
291
292
293
294  !*******************************************************************
295  ! Generate output: may be in concentration (ng/m3) or in mixing
296  ! ratio (ppt) or both
297  ! Output the position and the values alternated multiplied by
298  ! 1 or -1, first line is number of values, number of positions
299  ! For backward simulations, the unit is seconds, stored in grid_time
300  !*******************************************************************
301
302  ! Concentration output
303  !*********************
304  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
305
306  ! Wet deposition
307         sp_count_i=0
308         sp_count_r=0
309         sp_fact=-1.
310         sp_zer=.true.
311         if ((ldirect.eq.1).and.(WETDEP)) then
312         do jy=0,numygrid-1
313            do ix=0,numxgrid-1
314  !oncentraion greater zero
315              if (wetgrid(ix,jy).gt.smallnum) then
316                 if (sp_zer.eqv..true.) then ! first non zero value
317                    sp_count_i=sp_count_i+1
318                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
319                    sp_zer=.false.
320                    sp_fact=sp_fact*(-1.)
321                 endif
322                 sp_count_r=sp_count_r+1
323                 sparse_dump_r(sp_count_r)= &
324                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
325  !                sparse_dump_u(sp_count_r)=
326  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
327              else ! concentration is zero
328                  sp_zer=.true.
329              endif
330            end do
331         end do
332         else
333            sp_count_i=0
334            sp_count_r=0
335         endif
336         write(unitoutgrid) sp_count_i
337         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
338         write(unitoutgrid) sp_count_r
339         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
340  !       write(unitoutgrid) sp_count_u
341  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
342
343  ! Dry deposition
344         sp_count_i=0
345         sp_count_r=0
346         sp_fact=-1.
347         sp_zer=.true.
348         if ((ldirect.eq.1).and.(DRYDEP)) then
349          do jy=0,numygrid-1
350            do ix=0,numxgrid-1
351              if (drygrid(ix,jy).gt.smallnum) then
352                 if (sp_zer.eqv..true.) then ! first non zero value
353                    sp_count_i=sp_count_i+1
354                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
355                    sp_zer=.false.
356                    sp_fact=sp_fact*(-1.)
357                 endif
358                 sp_count_r=sp_count_r+1
359                 sparse_dump_r(sp_count_r)= &
360                      sp_fact* &
361                      1.e12*drygrid(ix,jy)/area(ix,jy)
362  !                sparse_dump_u(sp_count_r)=
363  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
364              else ! concentration is zero
365                  sp_zer=.true.
366              endif
367            end do
368          end do
369         else
370            sp_count_i=0
371            sp_count_r=0
372         endif
373         write(unitoutgrid) sp_count_i
374         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
375         write(unitoutgrid) sp_count_r
376         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
377  !       write(*,*) sp_count_u
378  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
379
380
381
382  ! Concentrations
383         sp_count_i=0
384         sp_count_r=0
385         sp_fact=-1.
386         sp_zer=.true.
387          do kz=1,numzgrid
388            do jy=0,numygrid-1
389              do ix=0,numxgrid-1
390                if (grid(ix,jy,kz).gt.smallnum) then
391                  if (sp_zer.eqv..true.) then ! first non zero value
392                    sp_count_i=sp_count_i+1
393                    sparse_dump_i(sp_count_i)= &
394                         ix+jy*numxgrid+kz*numxgrid*numygrid
395                    sp_zer=.false.
396                    sp_fact=sp_fact*(-1.)
397                   endif
398                   sp_count_r=sp_count_r+1
399                   sparse_dump_r(sp_count_r)= &
400                        sp_fact* &
401                        grid(ix,jy,kz)* &
402                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
403  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
404  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
405  !                sparse_dump_u(sp_count_r)=
406  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
407  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
408              else ! concentration is zero
409                  sp_zer=.true.
410              endif
411              end do
412            end do
413          end do
414         write(unitoutgrid) sp_count_i
415         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
416         write(unitoutgrid) sp_count_r
417         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
418  !       write(unitoutgrid) sp_count_u
419  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
420
421
422
423    endif !  concentration output
424
425  ! Mixing ratio output
426  !********************
427
428  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
429
430  ! Wet deposition
431         sp_count_i=0
432         sp_count_r=0
433         sp_fact=-1.
434         sp_zer=.true.
435         if ((ldirect.eq.1).and.(WETDEP)) then
436          do jy=0,numygrid-1
437            do ix=0,numxgrid-1
438                if (wetgrid(ix,jy).gt.smallnum) then
439                  if (sp_zer.eqv..true.) then ! first non zero value
440                    sp_count_i=sp_count_i+1
441                    sparse_dump_i(sp_count_i)= &
442                         ix+jy*numxgrid
443                    sp_zer=.false.
444                    sp_fact=sp_fact*(-1.)
445                 endif
446                 sp_count_r=sp_count_r+1
447                 sparse_dump_r(sp_count_r)= &
448                      sp_fact* &
449                      1.e12*wetgrid(ix,jy)/area(ix,jy)
450  !                sparse_dump_u(sp_count_r)=
451  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
452              else ! concentration is zero
453                  sp_zer=.true.
454              endif
455            end do
456          end do
457         else
458           sp_count_i=0
459           sp_count_r=0
460         endif
461         write(unitoutgridppt) sp_count_i
462         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
463         write(unitoutgridppt) sp_count_r
464         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
465  !       write(unitoutgridppt) sp_count_u
466  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
467
468
469  ! Dry deposition
470         sp_count_i=0
471         sp_count_r=0
472         sp_fact=-1.
473         sp_zer=.true.
474         if ((ldirect.eq.1).and.(DRYDEP)) then
475          do jy=0,numygrid-1
476            do ix=0,numxgrid-1
477                if (drygrid(ix,jy).gt.smallnum) then
478                  if (sp_zer.eqv..true.) then ! first non zero value
479                    sp_count_i=sp_count_i+1
480                    sparse_dump_i(sp_count_i)= &
481                         ix+jy*numxgrid
482                    sp_zer=.false.
483                    sp_fact=sp_fact*(-1)
484                 endif
485                 sp_count_r=sp_count_r+1
486                 sparse_dump_r(sp_count_r)= &
487                      sp_fact* &
488                      1.e12*drygrid(ix,jy)/area(ix,jy)
489  !                sparse_dump_u(sp_count_r)=
490  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
491              else ! concentration is zero
492                  sp_zer=.true.
493              endif
494            end do
495          end do
496         else
497           sp_count_i=0
498           sp_count_r=0
499         endif
500         write(unitoutgridppt) sp_count_i
501         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
502         write(unitoutgridppt) sp_count_r
503         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
504  !       write(unitoutgridppt) sp_count_u
505  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
506
507
508  ! Mixing ratios
509         sp_count_i=0
510         sp_count_r=0
511         sp_fact=-1.
512         sp_zer=.true.
513          do kz=1,numzgrid
514            do jy=0,numygrid-1
515              do ix=0,numxgrid-1
516                if (grid(ix,jy,kz).gt.smallnum) then
517                  if (sp_zer.eqv..true.) then ! first non zero value
518                    sp_count_i=sp_count_i+1
519                    sparse_dump_i(sp_count_i)= &
520                         ix+jy*numxgrid+kz*numxgrid*numygrid
521                    sp_zer=.false.
522                    sp_fact=sp_fact*(-1.)
523                 endif
524                 sp_count_r=sp_count_r+1
525                 sparse_dump_r(sp_count_r)= &
526                      sp_fact* &
527                      1.e12*grid(ix,jy,kz) &
528                      /volume(ix,jy,kz)/outnum* &
529                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
530  !                sparse_dump_u(sp_count_r)=
531  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
532  !+              outnum*weightair/weightmolar(ks)/
533  !+              densityoutgrid(ix,jy,kz)
534              else ! concentration is zero
535                  sp_zer=.true.
536              endif
537              end do
538            end do
539          end do
540         write(unitoutgridppt) sp_count_i
541         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
542         write(unitoutgridppt) sp_count_r
543         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
544  !       write(unitoutgridppt) sp_count_u
545  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
546
547      endif ! output for ppt
548
549  end do
550  end do
551
552    close(unitoutgridppt)
553    close(unitoutgrid)
554
555  end do
556
557  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
558  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
559       wetgridtotal
560  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
561       drygridtotal
562
563  ! Dump of receptor concentrations
564
565    if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
566      write(unitoutreceptppt) itime
567      do ks=1,nspec
568        write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
569             weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
570      end do
571    endif
572
573  ! Dump of receptor concentrations
574
575    if (numreceptor.gt.0) then
576      write(unitoutrecept) itime
577      do ks=1,nspec
578        write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
579             i=1,numreceptor)
580      end do
581    endif
582
583
584
585  ! Reinitialization of grid
586  !*************************
587
588  do ks=1,nspec
589  do kp=1,maxpointspec_act
590    do i=1,numreceptor
591      creceptor(i,ks)=0.
592    end do
593    do jy=0,numygrid-1
594      do ix=0,numxgrid-1
595        do l=1,nclassunc
596          do nage=1,nageclass
597            do kz=1,numzgrid
598              gridunc(ix,jy,kz,ks,kp,l,nage)=0.
599            end do
600          end do
601        end do
602      end do
603    end do
604  end do
605  end do
606
607end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG