source: trunk/MagicSoft/Mars/mhbase/MH3.cc@ 8238

Last change on this file since 8238 was 8082, checked in by tbretz, 18 years ago
*** empty log message ***
File size: 21.0 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 2002 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MH3
28//
29// With this histogram you can fill a histogram with up to three
30// variables from Mars parameter containers in an eventloop.
31//
32// In the constructor you can give up to three variables which should be
33// filled in the histogram. Dependend on the number of given variables
34// (data members) a TH1D, TH2D or TH3D is created.
35// Specify the data mamber like the following:
36// "MHillas.fLength"
37// Where MHillas is the name of the parameter container in the parameter
38// list and fLength is the name of the data member which should be filled
39// in the histogram. Assuming that your MHillas container has a different
40// name (MyHillas) the name to give would be:
41// "MyHillas.fLength"
42//
43// The axis binning is retrieved from the parameter list, too. Create a
44// MBinning with the name "Binning" plus the name of your MH3 container
45// plus the axis name ("X", "Y" or "Z") and add it to the parameter list.
46//
47// If you want to use a different unit for histogramming use SetScaleX,
48// SetScaleY and SetScaleZ.
49//
50//
51// Axis titles
52// ===========
53//
54// 1) If no other title is given the rule for this axis is used.
55// 2) If the MBinning used for this axis has a non-default Title
56// (MBinning::HasTitle) this title is used for the corresponding axis
57// 3) If the MH3 has a non-default title (MH3::SetTitle called)
58// this title is set as the histogram title. It can be used to overwrite
59// the axis titles. For more information see TH1::SetTitle, eg.
60// SetTitle("MyHist;x[mm];y[cm];Counts");
61//
62// For example:
63// MH3 myhist("MHillas.fLength");
64// myhist.SetName("MyHist");
65// myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
66// MBinning bins("BinningMyHistX");
67// bins.SetEdges(10, 0, 150);
68// plist.AddToList(&myhist);
69// plist.AddToList(&bins);
70//
71//
72// Class Version 2:
73// ----------------
74// - MDataChain *fData[3]; // Object from which the data is filled
75// + MData *fData[3]; // Object from which the data is filled
76//
77//
78/////////////////////////////////////////////////////////////////////////////
79#include "MH3.h"
80
81#include <ctype.h> // tolower
82#include <fstream>
83
84#include <TPad.h>
85#include <TStyle.h>
86#include <TCanvas.h>
87
88#include <TH2.h>
89#include <TH3.h>
90#include <TProfile.h>
91#include <TProfile2D.h>
92
93#include "MLog.h"
94#include "MLogManip.h"
95
96#include "MParList.h"
97#include "MBinning.h"
98#include "MDataPhrase.h"
99
100ClassImp(MH3);
101
102using namespace std;
103
104const TString MH3::gsDefName = "MH3";
105const TString MH3::gsDefTitle = "Container for a n-D Mars Histogram";
106
107// --------------------------------------------------------------------------
108//
109// Default constructor.
110//
111MH3::MH3(const unsigned int dim)
112 : fDimension(dim>3?3:dim), fHist(NULL)
113{
114 switch (fDimension)
115 {
116 case 1:
117 fHist = new TH1D;
118 fHist->SetYTitle("Counts");
119 break;
120 case 2:
121 fHist = new TH2D;
122 fHist->SetZTitle("Counts");
123 break;
124 case 3:
125 fHist = new TH3D;
126 break;
127 }
128
129 fData[0] = NULL;
130 fData[1] = NULL;
131 fData[2] = NULL;
132
133 fName = gsDefName;
134 fTitle = gsDefTitle;
135
136 if (fHist)
137 {
138 fHist->SetDirectory(NULL);
139 fHist->UseCurrentStyle();
140 }
141
142 fScale[0] = 1;
143 fScale[1] = 1;
144 fScale[2] = 1;
145}
146
147// --------------------------------------------------------------------------
148//
149// Creates an TH1D. memberx is filled into the X-bins. For a more detailed
150// description see the class description above.
151//
152MH3::MH3(const char *memberx)
153 : fDimension(1)
154{
155 fHist = new TH1D;
156
157 fData[0] = new MDataPhrase(memberx);
158 fData[1] = NULL;
159 fData[2] = NULL;
160
161 fName = gsDefName;
162 fTitle = gsDefTitle;
163
164 fHist->UseCurrentStyle();
165 fHist->SetDirectory(NULL);
166 fHist->SetYTitle("Counts");
167
168 fScale[0] = 1;
169 fScale[1] = 1;
170 fScale[2] = 1;
171}
172
173MH3::MH3(const TH1 &h1) : fDimension(1)
174{
175 if (h1.InheritsFrom(TH3::Class()))
176 fDimension = 3;
177 if (h1.InheritsFrom(TH2::Class()))
178 fDimension = 2;
179
180 fData[0] = NULL;
181 fData[1] = NULL;
182 fData[2] = NULL;
183
184 switch (fDimension)
185 {
186 case 3:
187 fData[2] = new MDataPhrase(h1.GetZaxis()->GetTitle());
188 case 2:
189 fData[1] = new MDataPhrase(h1.GetYaxis()->GetTitle());
190 case 1:
191 fData[0] = new MDataPhrase(h1.GetXaxis()->GetTitle());
192 }
193
194 fName = gsDefName;
195 fTitle = gsDefTitle;
196
197 fHist = (TH1*)h1.Clone();
198 fHist->SetDirectory(NULL);
199
200 fScale[0] = 1;
201 fScale[1] = 1;
202 fScale[2] = 1;
203}
204
205// --------------------------------------------------------------------------
206//
207// Creates an TH2D. memberx is filled into the X-bins. membery is filled
208// into the Y-bins. For a more detailed description see the class
209// description above.
210//
211MH3::MH3(const char *memberx, const char *membery)
212 : fDimension(2)
213{
214 fHist = new TH2D;
215
216 fData[0] = new MDataPhrase(memberx);
217 fData[1] = new MDataPhrase(membery);
218 fData[2] = NULL;
219
220 fName = gsDefName;
221 fTitle = gsDefTitle;
222
223 fHist->UseCurrentStyle();
224 fHist->SetDirectory(NULL);
225 fHist->SetZTitle("Counts");
226
227 fScale[0] = 1;
228 fScale[1] = 1;
229 fScale[2] = 1;
230}
231
232// --------------------------------------------------------------------------
233//
234// Creates an TH3D. memberx is filled into the X-bins. membery is filled
235// into the Y-bins. membery is filled into the Z-bins. For a more detailed
236// description see the class description above.
237//
238MH3::MH3(const char *memberx, const char *membery, const char *memberz)
239 : fDimension(3)
240{
241 fHist = new TH3D;
242
243 fData[0] = new MDataPhrase(memberx);
244 fData[1] = new MDataPhrase(membery);
245 fData[2] = new MDataPhrase(memberz);
246
247 fName = gsDefName;
248 fTitle = gsDefTitle;
249
250 fHist->UseCurrentStyle();
251 fHist->SetDirectory(NULL);
252
253 fScale[0] = 1;
254 fScale[1] = 1;
255 fScale[2] = 1;
256}
257
258// --------------------------------------------------------------------------
259//
260// Deletes the histogram
261//
262MH3::~MH3()
263{
264 delete fHist;
265
266 for (int i=0; i<3; i++)
267 if (fData[i])
268 delete fData[i];
269}
270
271// --------------------------------------------------------------------------
272//
273// Return the data members used by the data chain to be used in
274// MTask::AddBranchToList
275//
276TString MH3::GetDataMember() const
277{
278 TString str=fData[0]->GetDataMember();
279 if (fData[1])
280 {
281 str += ";";
282 str += fData[1]->GetDataMember();
283 }
284 if (fData[2])
285 {
286 str += ";";
287 str += fData[2]->GetDataMember();
288 }
289 return str;
290}
291
292// --------------------------------------------------------------------------
293//
294// Setup the Binning for the histograms automatically if the correct
295// instances of MBinning are found in the parameter list
296// For a more detailed description see class description above.
297//
298Bool_t MH3::SetupFill(const MParList *plist)
299{
300 // reset histogram (necessary if the same eventloop is run more than once)
301 fHist->Reset();
302
303 TString bname("Binning");
304 bname += fName;
305
306 MBinning *binsx = NULL;
307 MBinning *binsy = NULL;
308 MBinning *binsz = NULL;
309
310 switch (fDimension)
311 {
312 case 3:
313 binsz = (MBinning*)plist->FindObject(bname+"Z", "MBinning");
314 if (!binsz)
315 {
316 *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
317 return kFALSE;
318 }
319 if (fData[2] && !fData[2]->PreProcess(plist))
320 return kFALSE;
321 if (fData[2])
322 fHist->SetZTitle(fData[2]->GetTitle());
323 if (binsz->HasTitle())
324 fHist->SetZTitle(binsz->GetTitle());
325 if (binsz->IsLogarithmic())
326 fHist->SetBit(kIsLogz);
327 case 2:
328 binsy = (MBinning*)plist->FindObject(bname+"Y", "MBinning");
329 if (!binsy)
330 {
331 *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl;
332 return kFALSE;
333 }
334 if (fData[1] && !fData[1]->PreProcess(plist))
335 return kFALSE;
336 if (fData[1])
337 fHist->SetYTitle(fData[1]->GetTitle());
338 if (binsy->HasTitle())
339 fHist->SetYTitle(binsy->GetTitle());
340 if (binsy->IsLogarithmic())
341 fHist->SetBit(kIsLogy);
342 case 1:
343 binsx = (MBinning*)plist->FindObject(bname+"X", "MBinning");
344 if (!binsx)
345 {
346 if (fDimension==1)
347 binsx = (MBinning*)plist->FindObject(bname, "MBinning");
348
349 if (!binsx)
350 {
351 *fLog << err << dbginf << "Neither MBinning '" << bname << "X' nor '" << bname << "' found... aborting." << endl;
352 return kFALSE;
353 }
354 }
355 if (fData[0] && !fData[0]->PreProcess(plist))
356 return kFALSE;
357 if (fData[0]!=NULL)
358 fHist->SetXTitle(fData[0]->GetTitle());
359 if (binsx->HasTitle())
360 fHist->SetXTitle(binsx->GetTitle());
361 if (binsx->IsLogarithmic())
362 fHist->SetBit(kIsLogx);
363 }
364
365 TString title("Histogram for ");
366 title += fName;
367 title += Form(" (%dD)", fDimension);
368
369 fHist->SetName(fName);
370 fHist->SetTitle(fTitle==gsDefTitle ? title : fTitle);
371 fHist->SetDirectory(0);
372
373 switch (fDimension)
374 {
375 case 1:
376 SetBinning(fHist, binsx);
377 return kTRUE;
378 case 2:
379 SetBinning((TH2*)fHist, binsx, binsy);
380 return kTRUE;
381 case 3:
382 SetBinning((TH3*)fHist, binsx, binsy, binsz);
383 return kTRUE;
384 }
385
386 *fLog << err << "ERROR - MH3 has " << fDimension << " dimensions!" << endl;
387 return kFALSE;
388}
389
390// --------------------------------------------------------------------------
391//
392// Set the name of the histogram ant the MH3 container
393//
394void MH3::SetName(const char *name)
395{
396 if (fHist)
397 {
398 if (gPad)
399 {
400 const TString pfx(Form("%sProfX", fHist->GetName()));
401 const TString pfy(Form("%sProfY", fHist->GetName()));
402
403 TProfile *p = 0;
404 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
405 p->SetName(Form("%sProfX", name));
406 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
407 p->SetName(Form("%sProfY", name));
408 }
409
410 fHist->SetName(name);
411 fHist->SetDirectory(0);
412
413 }
414 MParContainer::SetName(name);
415}
416
417// --------------------------------------------------------------------------
418//
419// Set the title of the histogram ant the MH3 container
420//
421void MH3::SetTitle(const char *title)
422{
423 if (fHist)
424 fHist->SetTitle(title);
425 MParContainer::SetTitle(title);
426}
427
428// --------------------------------------------------------------------------
429//
430// Fills the one, two or three data members into our histogram
431//
432Bool_t MH3::Fill(const MParContainer *par, const Stat_t w)
433{
434 Double_t x=0;
435 Double_t y=0;
436 Double_t z=0;
437
438 switch (fDimension)
439 {
440 case 3:
441 z = fData[2]->GetValue()*fScale[2];
442 case 2:
443 y = fData[1]->GetValue()*fScale[1];
444 case 1:
445 x = fData[0]->GetValue()*fScale[0];
446 }
447
448 switch (fDimension)
449 {
450 case 3:
451 ((TH3*)fHist)->Fill(x, y, z, w);
452 return kTRUE;
453 case 2:
454 ((TH2*)fHist)->Fill(x, y, w);
455 return kTRUE;
456 case 1:
457 fHist->Fill(x, w);
458 return kTRUE;
459 }
460
461 return kFALSE;
462}
463/*
464// --------------------------------------------------------------------------
465//
466// Set the palette you wanna use:
467// - you could set the root "Pretty Palette Violet->Red" by
468// gStyle->SetPalette(1, 0), but in some cases this may look
469// confusing
470// - The maximum colors root allowes us to set by ourself
471// is 50 (idx: 51-100). This colors are set to a grayscaled
472// palette
473// - the number of contours must be two less than the number
474// of palette entries
475//
476void MHStarMap::PrepareDrawing() const
477{
478 const Int_t numg = 32; // number of gray scaled colors
479 const Int_t numw = 32; // number of white
480
481 Int_t palette[numg+numw];
482
483 //
484 // The first half of the colors are white.
485 // This is some kind of optical background supression
486 //
487 gROOT->GetColor(51)->SetRGB(1, 1, 1);
488
489 Int_t i;
490 for (i=0; i<numw; i++)
491 palette[i] = 51;
492
493 //
494 // now the (gray) scaled part is coming
495 //
496 for (;i<numw+numg; i++)
497 {
498 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
499
500 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
501 palette[i] = 52+i;
502 }
503
504 //
505 // Set the palette and the number of contour levels
506 //
507 gStyle->SetPalette(numg+numw, palette);
508 fStarMap->SetContour(numg+numw-2);
509}
510*/
511// --------------------------------------------------------------------------
512//
513// Setup a inversed deep blue sea palette for the fCenter histogram.
514//
515void MH3::SetColors() const
516{
517 // FIXME: This must be redone each time the canvas is repainted....
518 gStyle->SetPalette(51, NULL);
519 Int_t c[50];
520 for (int i=0; i<50; i++)
521 c[49-i] = gStyle->GetColorPalette(i);
522 gStyle->SetPalette(50, c);
523}
524
525void MH3::Paint(Option_t *o)
526{
527 TProfile *p=0;
528
529 const TString pfx(Form("%sProfX", fHist->GetName()));
530 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
531 {
532 Int_t col = p->GetLineColor();
533 p = ((TH2*)fHist)->ProfileX(pfx, -1, 9999, "s");
534 p->SetLineColor(col);
535 }
536
537 const TString pfy(Form("%sProfY", fHist->GetName()));
538 if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
539 {
540 Int_t col = p->GetLineColor();
541 p = ((TH2*)fHist)->ProfileY(pfy, -1, 9999, "s");
542 p->SetLineColor(col);
543 }
544
545 if (fHist->TestBit(kIsLogx) && fHist->GetEntries()>0) gPad->SetLogx();
546 if (fHist->TestBit(kIsLogy) && fHist->GetEntries()>0) gPad->SetLogy();
547 if (fHist->TestBit(kIsLogz) && fHist->GetEntries()>0) gPad->SetLogz();
548}
549
550// --------------------------------------------------------------------------
551//
552// Creates a new canvas and draws the histogram into it.
553//
554// Possible options are:
555// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
556// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
557// ONLY: Draw the profile histogram only (for 2D histograms only)
558// BLUE: Draw the profile in blue color instead of the histograms line color
559//
560// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
561// eg this is set when applying a logarithmic MBinning
562//
563// Be careful: The histogram belongs to this object and won't get deleted
564// together with the canvas.
565//
566void MH3::Draw(Option_t *opt)
567{
568 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
569 pad->SetBorderMode(0);
570 pad->SetGridx();
571 pad->SetGridy();
572
573 fHist->SetFillStyle(4000);
574
575 TString str(opt);
576 str.ToLower();
577
578 const Bool_t only = str.Contains("only") && fDimension==2;
579 const Bool_t same = str.Contains("same") && fDimension==2;
580 const Bool_t blue = str.Contains("blue") && fDimension==2;
581 const Bool_t profx = str.Contains("profx") && fDimension==2;
582 const Bool_t profy = str.Contains("profy") && fDimension==2;
583
584 str.ReplaceAll("only", "");
585 str.ReplaceAll("blue", "");
586 str.ReplaceAll("profx", "");
587 str.ReplaceAll("profy", "");
588 str.ReplaceAll(" ", "");
589
590 // FIXME: We may have to remove all our own options from str!
591 if (!only)
592 fHist->Draw(str);
593
594 AppendPad();
595
596 TProfile *p=0;
597 if (profx)
598 {
599 const TString pfx(Form("%sProfX", fHist->GetName()));
600
601 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
602 *fLog << warn << "TProfile " << pfx << " already in pad." << endl;
603
604 p = ((TH2*)fHist)->ProfileX(pfx, -1, 9999, "s");
605 p->UseCurrentStyle();
606 p->SetLineColor(blue ? kBlue : fHist->GetLineColor());
607 p->SetBit(kCanDelete);
608 p->SetDirectory(NULL);
609 p->SetXTitle(fHist->GetXaxis()->GetTitle());
610 p->SetYTitle(fHist->GetYaxis()->GetTitle());
611 p->Draw(only&&!same?"":"same");
612 }
613 if (profy)
614 {
615 const TString pfy(Form("%sProfY", fHist->GetName()));
616
617 if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
618 *fLog << warn << "TProfile " << pfy << " already in pad." << endl;
619
620 p = ((TH2*)fHist)->ProfileY(pfy, -1, 9999, "s");
621 p->UseCurrentStyle();
622 p->SetLineColor(blue ? kBlue : fHist->GetLineColor());
623 p->SetBit(kCanDelete);
624 p->SetDirectory(NULL);
625 p->SetYTitle(fHist->GetXaxis()->GetTitle());
626 p->SetXTitle(fHist->GetYaxis()->GetTitle());
627 p->Draw(only&&!same?"":"same");
628 }
629
630 //AppendPad("log");
631}
632
633// --------------------------------------------------------------------------
634//
635// Implementation of SavePrimitive. Used to write the call to a constructor
636// to a macro. In the original root implementation it is used to write
637// gui elements to a macro-file.
638//
639void MH3::StreamPrimitive(ostream &out) const
640{
641 TString name = GetUniqueName();
642
643 out << " MH3 " << name << "(\"";
644 out << fData[0]->GetRule() << "\"";
645 if (fDimension>1)
646 out << ", \"" << fData[1]->GetRule() << "\"";
647 if (fDimension>2)
648 out << ", \"" << fData[2]->GetRule() << "\"";
649
650 out << ");" << endl;
651
652 if (fName!=gsDefName)
653 out << " " << name << ".SetName(\"" << fName << "\");" << endl;
654
655 if (fTitle!=gsDefTitle)
656 out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
657
658 switch (fDimension)
659 {
660 case 3:
661 if (fScale[2]!=1)
662 out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
663 case 2:
664 if (fScale[1]!=1)
665 out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
666 case 1:
667 if (fScale[0]!=1)
668 out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
669 }
670}
671
672// --------------------------------------------------------------------------
673//
674// Used to rebuild a MH3 object of the same type (data members,
675// dimension, ...)
676//
677MParContainer *MH3::New() const
678{
679 MH3 *h = NULL;
680
681 if (fData[0] == NULL)
682 h=new MH3(fDimension);
683 else
684 switch (fDimension)
685 {
686 case 1:
687 h=new MH3(fData[0]->GetRule());
688 break;
689 case 2:
690 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule());
691 break;
692 case 3:
693 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule());
694 break;
695 }
696 switch (fDimension)
697 {
698 case 3:
699 h->SetScaleZ(fScale[2]);
700 case 2:
701 h->SetScaleY(fScale[1]);
702 case 1:
703 h->SetScaleX(fScale[0]);
704 }
705 return h;
706}
707
708TString MH3::GetRule(const Char_t axis) const
709{
710 switch (tolower(axis))
711 {
712 case 'x':
713 return fData[0] ? fData[0]->GetRule() : TString("");
714 case 'y':
715 return fData[1] ? fData[1]->GetRule() : TString("");
716 case 'z':
717 return fData[2] ? fData[2]->GetRule() : TString("");
718 default:
719 return "<n/a>";
720 }
721}
722
723// --------------------------------------------------------------------------
724//
725// Returns the total number of bins in a histogram (excluding under- and
726// overflow bins)
727//
728Int_t MH3::GetNbins() const
729{
730 Int_t num = 1;
731
732 switch (fDimension)
733 {
734 case 3:
735 num *= fHist->GetNbinsZ()+2;
736 case 2:
737 num *= fHist->GetNbinsY()+2;
738 case 1:
739 num *= fHist->GetNbinsX()+2;
740 }
741
742 return num;
743}
744
745// --------------------------------------------------------------------------
746//
747// Returns the total number of bins in a histogram (excluding under- and
748// overflow bins) Return -1 if bin is underflow or overflow
749//
750Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
751{
752 const TAxis &axex = *fHist->GetXaxis();
753 const TAxis &axey = *fHist->GetYaxis();
754 const TAxis &axez = *fHist->GetZaxis();
755
756 Int_t binz = 0;
757 Int_t biny = 0;
758 Int_t binx = 0;
759
760 switch (fDimension)
761 {
762 case 3:
763 binz = axez.FindFixBin(z);
764 if (binz>axez.GetLast() || binz<axez.GetFirst())
765 return -1;
766 case 2:
767 biny = axey.FindFixBin(y);
768 if (biny>axey.GetLast() || biny<axey.GetFirst())
769 return -1;
770 case 1:
771 binx = axex.FindFixBin(x);
772 if (binx<axex.GetFirst() || binx>axex.GetLast())
773 return -1;
774 }
775
776 const Int_t nx = fHist->GetNbinsX()+2;
777 const Int_t ny = fHist->GetNbinsY()+2;
778
779 return binx + nx*(biny +ny*binz);
780}
Note: See TracBrowser for help on using the repository browser.