KaliVeda  1.12/06
Heavy-Ion Analysis Toolkit
KVSimDir.cpp
Go to the documentation of this file.
1 //Created by KVClassFactory on Tue Jul 17 08:46:41 2012
2 //Author: John Frankland,,,
3 
4 #include "KVSimDir.h"
5 #include "KVSimFile.h"
6 #include "TSystemDirectory.h"
7 #include "TSystem.h"
8 #include "KVString.h"
9 #include "TTree.h"
10 #include "TKey.h"
11 #include "TFile.h"
12 #include "TBranchElement.h"
13 #include <iostream>
14 using namespace std;
15 
17 
18 
19 
20 
21 
25 {
26  // Default constructor
27  init();
28 }
29 
30 
31 
34 
35 KVSimDir::KVSimDir(const Char_t* name, const Char_t* path)
36  : KVBase(name, path)
37 {
38  // Ctor with name and directory to analyse
39  init();
40  // make sure we have absolute path
41  KVString _path(path);
42  gSystem->ExpandPathName(_path);
43  if (!gSystem->IsAbsoluteFileName(_path)) {
44  if (_path.CountChar('.') == 1) {
45  _path.ReplaceAll(".", "$(PWD)");
46  gSystem->ExpandPathName(_path);
47  }
48  }
49  SetTitle(_path);
50 }
51 
52 
53 
56 
58 {
59  // Default initialisations
60  fSimData.SetName("Simulated Data");
61  fFiltData.SetName("Filtered Simulated Data");
62 }
63 
64 
65 
66 
73 
75 {
76  // Copy constructor
77  // This ctor is used to make a copy of an existing object (for example
78  // when a method returns an object), and it is always a good idea to
79  // implement it.
80  // If your class allocates memory in its constructor(s) then it is ESSENTIAL :-)
81 
82  obj.Copy(*this);
83 }
84 
85 
86 
89 
91 {
92  // Destructor
93 }
94 
95 
96 
97 
105 
106 void KVSimDir::Copy(TObject& obj) const
107 {
108  // This method copies the current state of 'this' object into 'obj'
109  // You should add here any member variables, for example:
110  // (supposing a member variable KVSimDir::fToto)
111  // CastedObj.fToto = fToto;
112  // or
113  // CastedObj.SetToto( GetToto() );
114 
115  KVBase::Copy(obj);
116  //KVSimDir& CastedObj = (KVSimDir&)obj;
117 }
118 
119 
120 
121 
125 
127 {
128  // Read contents of directory given to ctor.
129  // Each ROOT file will be analysed by AnalyseFile().
130 
131  Info("AnalyseDirectory", "Analysing %s...", GetDirectory());
132  fSimData.Clear();
133  fFiltData.Clear();
134  //loop over files in current directory
135  TSystemDirectory thisDir(".", GetDirectory());
136  unique_ptr<TList> fileList(thisDir.GetListOfFiles());
137  TIter nextFile(fileList.get());
138  TSystemFile* aFile = 0;
139  while ((aFile = (TSystemFile*)nextFile())) {
140 
141  if (aFile->IsDirectory()) continue;
142 
143  KVString fileName = aFile->GetName();
144  if (!fileName.EndsWith(".root") && !fileName.Contains(".root.")) continue; /* skip non-ROOT files */
145 
146  AnalyseFile(fileName);
147 
148  }
149 }
150 
151 
152 
153 
173 
174 void KVSimDir::AnalyseFile(const Char_t* filename)
175 {
176  // Analyse ROOT file given as argument.
177  // If there is a TTree in the file, then we look at all of its branches until we find one
178  // containing objects which derive from KVEvent:
179  //
180  // -- if they inherit from KVSimEvent, we add the file to the list of simulated data:
181  // * a KVSimFile is created. The title of the TTree where data were found will
182  // be used as 'Information' on the nature of the simulation.
183  // -- if they inherit from KVReconstructedEvent, we add the file to the list of filtered data.
184  // * a KVSimFile is created. Informations on the filtered data are extracted from
185  // TNamed objects in the file with names 'Dataset', 'System', 'Run', 'Geometry'
186  // (type of geometry used, 'ROOT' or 'KV'), 'Origin' (i.e. the name of the simulation
187  // file which was filtered), 'Filter' (type of filter: Geo, GeoThresh or Full).
188  // These objects are automatically created when data is filtered using KVEventFiltering.
189  //
190  // Analysis of the file stops after the first TTree with a branch satisfying one of the
191  // two criteria is found (it is assumed that in each file there is only one TTree containing
192  // either simulated or filtered data).
193 
194  //Info("AnalyseFile", "Analysing file %s...", filename);
195  TString fullpath;
196  AssignAndDelete(fullpath, gSystem->ConcatFileName(GetDirectory(), filename));
197  unique_ptr<TFile> file(TFile::Open(fullpath));
198  if (!file.get() || file->IsZombie()) return;
199  // look for TTrees in file
200  TIter next(file->GetListOfKeys());
201  TKey* key;
202  while ((key = (TKey*)next())) {
203  TString cn = key->GetClassName();
204  if (cn == "TTree") {
205  // look for branch with KVEvent objects
206  TTree* tree = (TTree*)file->Get(key->GetName());
207  TSeqCollection* branches = tree->GetListOfBranches();
208  TIter nextB(branches);
209  TBranchElement* branch;
210  while ((branch = (TBranchElement*)nextB())) {
211  TString branch_classname = branch->GetClassName();
212  TClass* branch_class = TClass::GetClass(branch_classname, kFALSE, kTRUE);
213  if (branch_class && branch_class->InheritsFrom("KVEvent")) {
214  if (branch_class->InheritsFrom("KVReconstructedEvent")) {
215  // filtered data. there must be TNamed called 'Dataset', 'System', & 'Run' in the file.
216  unique_ptr<TNamed> ds((TNamed*)file->Get("Dataset"));
217  unique_ptr<TNamed> orig((TNamed*)file->Get("Origin"));
218  unique_ptr<TNamed> sys((TNamed*)file->Get("System"));
219  unique_ptr<TNamed> r((TNamed*)file->Get("Run"));
220  unique_ptr<TNamed> g((TNamed*)file->Get("Geometry"));
221  unique_ptr<TNamed> f((TNamed*)file->Get("Filter"));
222  TString dataset;
223  if (ds.get()) dataset = ds->GetTitle();
224  TString system;
225  if (sys.get()) system = sys->GetTitle();
226  TString run;
227  if (r.get()) run = r->GetTitle();
228  TString origin;
229  if (orig.get()) origin = orig->GetTitle();
231  if (g.get()) geometry = g->GetTitle();
232  TString filter;
233  if (f.get()) filter = f->GetTitle();
234  Int_t run_number = run.Atoi();
235  KVSimFile* fff = new KVSimFile(this, filename, tree->GetTitle(), tree->GetEntries(), tree->GetName(), branch->GetName(),
236  dataset, system, run_number, geometry, origin, filter);
237  fFiltData.Add(fff);
238  unique_ptr<TNamed> gem((TNamed*)file->Get("Gemini++"));
239  if (gem.get()) {
240  if (!strcmp(gem->GetTitle(), "yes")) fff->SetGemini();
241  unique_ptr<TNamed> gemdec((TNamed*)file->Get("GemDecayPerEvent"));
242  if (gemdec.get()) {
243  TString gemdecperev = gemdec->GetTitle();
244  fff->SetGemDecayPerEvent(gemdecperev.Atoi());
245  }
246  }
247  return;
248  }
249  else {
250  fSimData.Add(new KVSimFile(this, filename, tree->GetTitle(), tree->GetEntries(), tree->GetName(), branch->GetName()));
251  return;
252  }
253  }
254  }
255  }
256  }
257 }
258 
259 
260 
262 
264 {
265  fSimData.Add(f);
266 }
267 
268 
269 
271 
273 {
274  fFiltData.Add(f);
275 }
276 
277 
278 
279 
281 
282 void KVSimDir::ls(Option_t*) const
283 {
285  cout << "SIMULATION SET: " << GetName() << endl;
286  cout << "DIRECTORY: " << GetDirectory() << endl;
287  cout << "CONTENTS:" << endl;
288  cout << "--simulated data:" << endl;
289  fSimData.ls();
290  cout << "--filtered data:" << endl;
291  fFiltData.ls();
292 }
293 
294 
int Int_t
void AssignAndDelete(TString &target, char *tobedeleted)
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
ROOT::R::TRInterface & r
#define f(i)
char Char_t
const Bool_t kFALSE
const Bool_t kTRUE
const char Option_t
R__EXTERN TSystem * gSystem
Base class for KaliVeda framework.
Definition: KVBase.h:135
virtual void Copy(TObject &) const
Make a copy of this object.
Definition: KVBase.cpp:397
virtual void Clear(Option_t *option="")
virtual void Add(TObject *obj)
Handle directory containing simulated and/or filtered simulated data ,.
Definition: KVSimDir.h:43
void init()
Default initialisations.
Definition: KVSimDir.cpp:57
KVSimDir()
Default constructor.
Definition: KVSimDir.cpp:24
void ls(Option_t *opt="") const
Definition: KVSimDir.cpp:282
virtual void AnalyseFile(const Char_t *)
Definition: KVSimDir.cpp:174
virtual ~KVSimDir()
Destructor.
Definition: KVSimDir.cpp:90
KVList fSimData
list of simulated data files
Definition: KVSimDir.h:45
KVList fFiltData
list of filtered simulated data files
Definition: KVSimDir.h:46
void AddSimData(KVSimFile *)
Definition: KVSimDir.cpp:263
virtual void AnalyseDirectory()
Definition: KVSimDir.cpp:126
void AddFiltData(KVSimFile *)
Definition: KVSimDir.cpp:272
virtual const Char_t * GetDirectory() const
Definition: KVSimDir.h:60
void Copy(TObject &) const
Definition: KVSimDir.cpp:106
Handle file containing simulated and/or filtered simulated data ,.
Definition: KVSimFile.h:18
void SetGemini(Bool_t yes=kTRUE)
Definition: KVSimFile.h:52
void SetGemDecayPerEvent(Int_t n)
Definition: KVSimFile.h:56
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
virtual const char * GetClassName() const
static TClass * GetClass(Bool_t load=kTRUE, Bool_t silent=kFALSE)
Bool_t InheritsFrom(const char *cl) const
virtual void ls(Option_t *option="") const
void SetName(const char *name)
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
virtual const char * GetClassName() const
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual const char * GetTitle() const
virtual void Info(const char *method, const char *msgfmt,...) const
static void IndentLevel()
Int_t Atoi() const
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Int_t CountChar(Int_t c) const
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
TString & ReplaceAll(const char *s1, const char *s2)
virtual TList * GetListOfFiles() const
virtual Bool_t IsDirectory(const char *dir=0) const
virtual char * ConcatFileName(const char *dir, const char *name)
virtual Bool_t IsAbsoluteFileName(const char *dir)
virtual char * ExpandPathName(const char *path)
const long double g
masses
Definition: KVUnits.h:72