StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StHistCollectorMaker.cxx
1 //*-- Author : Valeri Fine (fine@bnl.gov)
2 //
3 // $Id: StHistCollectorMaker.cxx,v 2.7 2007/05/29 20:46:19 fine Exp $
4 //
6 // //
7 // StHistCollectorMaker is to collect the histBranch staff from the several //
8 // files //
9 // //
10 // use $STAR/StRoot/macro/analysis/doHists.C macro to see hot it works //
11 // //
12 // This maker collects non-empty histograms from histBranch //
13 // calling GetDataSet("hist") //
14 // and accumulates it under ".const/Merged" sub-dataset. //
15 // The "Merged" dataset preserves the dataset structure of the original //
16 // histBrach. Namely, it contains the one level of heirarchy with //
17 // the maker names and TList object contaoins the list of histogram //
18 // producing by the maker with the several reconstruction sessions //
19 // //
20 // See: $STAR/StRoot/macros/analysis/doHists.C macro as example. //
21 // doHist.C macro is derived from doEvents.C. It preserves doEvents.C
22 // user interface.
23 // //
24 // Submit any problem with this code via begin_html <A HREF="http://www.star.bnl.gov/STARAFS/comp/sofi/bugs/send-pr.html"><B><I>"STAR Problem Report Form"</I></B></A> end_html //
25 // //
27 
28 #include "StHistCollectorMaker.h"
29 #include "TDataSetIter.h"
30 #include "TH1.h"
31 
32 ClassImp(StHistCollectorMaker)
33 
34 //_____________________________________________________________________________
35 StHistCollectorMaker::StHistCollectorMaker(const char *name):StMaker(name){
36  // StHistCollectorMaker constructor
37  // const char *name - the name of the maker
38  fMergedSet = 0;
39 }
40 //_____________________________________________________________________________
41 StHistCollectorMaker::~StHistCollectorMaker(){ }
42 //_____________________________________________________________________________
43 Int_t StHistCollectorMaker::Init(){
44 // Create the dataset to merger the coming histograms
45  if (fMergedSet) delete fMergedSet;
46  fMergedSet = 0;
47  return StMaker::Init();
48 }
49 //_____________________________________________________________________________
51  // Make - this methoid is called in loop for each event
52  static const char *datasetName = "Merged";
53  if (!fMergedSet) {
54  // Qurey one's "Mergeg" dataset (assuming I/O maker)
55  // to take over
56  fMergedSet = GetDataSet(datasetName);
57  if (!fMergedSet) fMergedSet = new TDataSet(datasetName);
58  fMergedSet->Shunt(0);
59  AddConst(fMergedSet);
60  }
61  AddHists();
62  return kStOK;
63 }
64 
65 //_____________________________________________________________________________
66 TDataSet *StHistCollectorMaker::AddHists()
67 {
68  // Update dataset fMergedSet with "histBranch" staff
69  TDataSet *rec = fMergedSet;
70  TDataSet *histBranch = GetDataSet("hist");
71  if (!histBranch) {
72  LOG_INFO << "No hist" << endm;
73  return 0;}
74  TDataSetIter nextDonor(histBranch);
75  TDataSet *donor = 0;
76  while((donor = nextDonor())) {
77  Bool_t found = kFALSE;
78  TDataSetIter nextOld(rec);
79  const Char_t *newname = donor->GetName();
80  TDataSet *oldset = 0;
81  while ( (oldset = nextOld()) && !found) {
82  // if the "new" set does contain the dataset
83  // with the same name as ours update it too
84  if (oldset->IsThisDir(newname)) {
85  UpdateHists((TObjectSet *)oldset,(TObjectSet *)donor);
86  found = kTRUE;
87  }
88  }
89  // If the new "set" contains some new dataset with brand-new name
90  // move it into the our dataset and remove it from its old location
91  if (!found) {
92  // Take in account the non-empty list only
93  TList *newList = (TList *)donor->GetObject();
94  if (newList && newList->GetSize() > 0) {
95  // remove the histogram from its directories
96  TIter nextDonor(newList);
97  TH1 *donorHist = 0;
98  Int_t count = 0;
99  while ( (donorHist = (TH1 *)nextDonor()) ) {
100  if (donorHist->GetEntries() > 0) {
101  donorHist->SetDirectory(0);
102  count++; // count non-empty histograms
103  } else {
104  newList->Remove(donorHist);
105  }
106  }
107  if (count) donor->Shunt(rec); // accept this if it contains one non-empty hist at least
108  }
109  }
110  }
111  return 0;
112 }
113 //_____________________________________________________________________________
114 void StHistCollectorMaker::UpdateHists(TObjectSet *oldSet,TObjectSet *newSet)
115 {
116  // Merge two sets of the histograms
117  if (oldSet && newSet) {
118  // Get the list of the new histograms
119  TList *newList = (TList *)newSet->GetObject();
120  TList *oldList = (TList *)oldSet->GetObject();
121  TIter nextDonor(newList);
122  TIter nextOld(oldList);
123 
124  TH1 *donor = 0;
125  while ( (donor = (TH1 *)nextDonor()) ) {
126  Bool_t found = kFALSE;
127  if (donor->GetEntries() > 0 ) {
128  const Char_t *newname = donor->GetName();
129  TH1 *oldHist = 0;
130  while ( (oldHist = (TH1 *)nextOld()) && !found) {
131  // if the "new" set does contain the dataset
132  // with the same name as ours update it too
133  if (!strcmp(oldHist->GetName(),newname)) {
134  found = kTRUE;
135  oldHist->Add(donor); // merge histograms
136  }
137  }; nextOld.Reset(); // old histograms has been looked up
138 
139  // If the new "set" contains some new dataset with brand-new name
140  // move it into the our dataset and remove it from its old location
141  if (!found) { // move histogram to this
142  oldList->Add(donor);
143  newList->Remove(donor);
144  donor->SetDirectory(0);
145  }
146  }
147  }
148  }
149 }
150 
151 // $Log: StHistCollectorMaker.cxx,v $
152 // Revision 2.7 2007/05/29 20:46:19 fine
153 // Introduce logger-based output
154 //
155 // Revision 2.6 2007/04/28 20:36:14 perev
156 // Redundant StChain.h removed
157 //
158 // Revision 2.5 2000/12/02 01:03:10 fine
159 // Recursive collection of histogram introduced
160 //
161 // Revision 2.4 2000/12/01 02:16:42 fine
162 // the performance issue resolved
163 //
164 // Revision 2.3 2000/11/30 20:13:27 fine
165 // Get rid of the empty histograms (thanks Fisyak)
166 //
167 // Revision 2.2 2000/11/30 19:37:26 fine
168 // Reference to doHists.C macro has been added
169 //
170 // Revision 2.1 2000/11/30 19:35:14 fine
171 // New analysis utility to collect all histogram from all histBranh production branches
172 //
173 
virtual Bool_t IsThisDir(const char *dirname, int len=-1, int ignorecase=0) const
Definition: TDataSet.cxx:555
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:40
virtual void Shunt(TDataSet *newParent=0)
Definition: TDataSet.cxx:810
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TObjectSet.h:56