Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TSelEventGen.cxx
Go to the documentation of this file.
1 // @(#)root/proof:$Id$
2 // Author: Sangsu Ryu 22/06/2010
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2005, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 /** \class TSelEventGen
13 \ingroup proofbench
14 
15 Selector for event file generation.
16 List of files to be generated for each node is provided by client.
17 And list of files generated is sent back.
18 Existing files are reused if not forced to be regenerated.
19 
20 */
21 
22 #define TSelEventGen_cxx
23 
24 #include "TSelEventGen.h"
25 #include "TParameter.h"
26 #include "TProofNodeInfo.h"
27 #include "TProofBenchTypes.h"
28 #include "TProof.h"
29 #include "TMap.h"
30 #include "TDSet.h"
31 #include "TEnv.h"
32 #include "TFileInfo.h"
33 #include "TFile.h"
34 #include "TSortedList.h"
35 #include "TRandom.h"
36 #include "Event.h"
37 #include "TProofServ.h"
38 #include "TMacro.h"
39 
40 ClassImp(TSelEventGen);
41 
42 ////////////////////////////////////////////////////////////////////////////////
43 /// Constructor
44 
45 TSelEventGen::TSelEventGen()
46  : fBaseDir(""), fNEvents(100000), fNTracks(100), fNTracksMax(-1),
47  fRegenerate(kFALSE), fTotalGen(0), fFilesGenerated(0),
48  fGenerateFun(0), fChain(0)
49 {
50  if (gProofServ){
51  fBaseDir=gProofServ->GetDataDir();
52  // Two directories up
53  fBaseDir.Remove(fBaseDir.Last('/'));
54  fBaseDir.Remove(fBaseDir.Last('/'));
55  }
56  else{
57  fBaseDir="";
58  }
59 }
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 /// The Begin() function is called at the start of the query.
63 /// When running with PROOF Begin() is only called on the client.
64 /// The tree argument is deprecated (on PROOF 0 is passed).
65 
66 void TSelEventGen::Begin(TTree *)
67 {
68  TString option = GetOption();
69  // Determine the test type
70  TMap *filemap = dynamic_cast<TMap *>
71  (fInput->FindObject("PROOF_FilesToProcess"));
72  if (filemap) {
73  //Info("Begin", "dumping the file map:");
74  //filemap->Print();
75  } else {
76  if (fInput->FindObject("PROOF_FilesToProcess")) {
77  Error("Begin", "object 'PROOF_FilesToProcess' found but not a map"
78  " (%s)", fInput->FindObject("PROOF_FilesToProcess")->ClassName());
79  } else {
80  Error("Begin", "object 'PROOF_FilesToProcess' not found");
81  }
82  }
83 }
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 /// The SlaveBegin() function is called after the Begin() function.
87 /// When running with PROOF SlaveBegin() is called on each slave server.
88 /// The tree argument is deprecated (on PROOF 0 is passed).
89 
90 void TSelEventGen::SlaveBegin(TTree *tree)
91 {
92  Init(tree);
93 
94  TString option = GetOption();
95 
96  //get parameters
97 
98  Bool_t found_basedir=kFALSE;
99  Bool_t found_nevents=kFALSE;
100  Bool_t found_ntracks=kFALSE;
101  Bool_t found_ntrkmax=kFALSE;
102  Bool_t found_regenerate=kFALSE;
103 
104  TIter nxt(fInput);
105  TString sinput;
106  TObject *obj;
107 
108  while ((obj = nxt())){
109 
110  sinput=obj->GetName();
111  //Info("SlaveBegin", "Input list: %s", sinput.Data());
112 
113  if (sinput.Contains("PROOF_BenchmarkBaseDir")){
114  TNamed* a = dynamic_cast<TNamed*>(obj);
115  if (a){
116  TString bdir = a->GetTitle();
117  if (!bdir.IsNull()){
118  TUrl u(bdir, kTRUE);
119  Bool_t isLocal = !strcmp(u.GetProtocol(), "file") ? kTRUE : kFALSE;
120  if (isLocal && !gSystem->IsAbsoluteFileName(u.GetFile()))
121  u.SetFile(TString::Format("%s/%s", fBaseDir.Data(), u.GetFile()));
122  if (isLocal) {
123  if ((gSystem->AccessPathName(u.GetFile()) &&
124  gSystem->mkdir(u.GetFile(), kTRUE) == 0) ||
125  !gSystem->AccessPathName(u.GetFile(), kWritePermission)) {
126  // Directory is writable
127  fBaseDir = u.GetFile();
128  Info("SlaveBegin", "using base directory \"%s\"", fBaseDir.Data());
129  } else{
130  // Directory is not writable or not available, use default directory
131  Warning("SlaveBegin", "\"%s\" directory is not writable or not existing,"
132  " using default directory: %s",
133  bdir.Data(), fBaseDir.Data());
134  }
135  } else {
136  // We assume the user knows what it does
137  fBaseDir = bdir;
138  Info("SlaveBegin", "using non local base directory \"%s\"", fBaseDir.Data());
139  }
140  } else{
141  Info("SlaveBegin", "using default directory: %s",
142  fBaseDir.Data());
143  }
144  found_basedir=kTRUE;
145  }
146  else{
147  Error("SlaveBegin", "PROOF_BenchmarkBaseDir not type TNamed");
148  }
149  continue;
150  }
151  if (sinput.Contains("PROOF_BenchmarkNEvents")){
152  TParameter<Long64_t>* a=dynamic_cast<TParameter<Long64_t>*>(obj);
153  if (a){
154  fNEvents= a->GetVal();
155  found_nevents=kTRUE;
156  }
157  else{
158  Error("SlaveBegin", "PROOF_BenchmarkEvents not type TParameter"
159  "<Long64_t>*");
160  }
161  continue;
162  }
163  if (sinput.Contains("PROOF_BenchmarkNTracks")){
164  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
165  if (a){
166  fNTracks=a->GetVal();
167  found_ntracks=kTRUE;
168  }
169  else{
170  Error("SlaveBegin", "PROOF_BenchmarkNTracks not type TParameter"
171  "<Int_t>*");
172  }
173  continue;
174  }
175  if (sinput.Contains("PROOF_BenchmarkNTracksMax")){
176  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
177  if (a){
178  fNTracksMax=a->GetVal();
179  found_ntrkmax=kTRUE;
180  }
181  else{
182  Error("SlaveBegin", "PROOF_BenchmarkNTracksMax not type TParameter"
183  "<Int_t>*");
184  }
185  continue;
186  }
187  if (sinput.Contains("PROOF_BenchmarkRegenerate")){
188  TParameter<Int_t>* a=dynamic_cast<TParameter<Int_t>*>(obj);
189  if (a){
190  fRegenerate=a->GetVal();
191  found_regenerate=kTRUE;
192  }
193  else{
194  Error("SlaveBegin", "PROOF_BenchmarkRegenerate not type TParameter"
195  "<Int_t>*");
196  }
197  continue;
198  }
199  if (sinput.Contains("PROOF_GenerateFun")){
200  TNamed *a = dynamic_cast<TNamed*>(obj);
201  if (!(fGenerateFun = dynamic_cast<TMacro *>(fInput->FindObject(a->GetTitle())))) {
202  Error("SlaveBegin", "PROOF_GenerateFun requires the TMacro object in the input list");
203  }
204  continue;
205  }
206  }
207 
208  if (!found_basedir){
209  Warning("SlaveBegin", "PROOF_BenchmarkBaseDir not found; using default:"
210  " %s", fBaseDir.Data());
211  }
212  if (!found_nevents){
213  Warning("SlaveBegin", "PROOF_BenchmarkNEvents not found; using default:"
214  " %lld", fNEvents);
215  }
216  if (!found_ntracks){
217  Warning("SlaveBegin", "PROOF_BenchmarkNTracks not found; using default:"
218  " %d", fNTracks);
219  }
220  if (!found_ntrkmax){
221  Warning("SlaveBegin", "PROOF_BenchmarkNTracksMax not found; using default:"
222  " %d", fNTracksMax);
223  } else if (fNTracksMax <= fNTracks) {
224  Warning("SlaveBegin", "PROOF_BenchmarkNTracksMax must be larger then"
225  " fNTracks=%d ; ignoring", fNTracks);
226  fNTracksMax = -1;
227  found_ntrkmax = kFALSE;
228  }
229  if (!found_regenerate){
230  Warning("SlaveBegin", "PROOF_BenchmarkRegenerate not found; using"
231  " default: %d", fRegenerate);
232  }
233 
234  fFilesGenerated = new TList();
235 
236  TString hostname(TUrl(gSystem->HostName()).GetHostFQDN());
237  TString thisordinal = gProofServ ? gProofServ->GetOrdinal() : "n.d";
238  TString sfilegenerated =
239  TString::Format("PROOF_FilesGenerated_%s_%s", hostname.Data(), thisordinal.Data());
240  fFilesGenerated->SetName(sfilegenerated);
241  fFilesGenerated->SetOwner(kTRUE);
242 }
243 
244 ////////////////////////////////////////////////////////////////////////////////
245 ///Generate files for IO-bound run
246 ///Input parameters
247 /// filename: The name of the file to be generated
248 /// sizenevents: Either the number of events to generate when
249 /// filetype==kPBFileBenchmark
250 /// or the size of the file to generate when
251 /// filetype==kPBFileCleanup
252 ///Returns
253 /// Either Number of entries in the file when
254 /// filetype==kPBFileBenchmark
255 /// or bytes written when filetype==kPBFileCleanup
256 ///return 0 in case error
257 
258 Long64_t TSelEventGen::GenerateFiles(const char *filename, Long64_t sizenevents)
259 {
260  Long64_t nentries=0;
261  TDirectory* savedir = gDirectory;
262  //printf("current dir=%s\n", gDirectory->GetPath());
263 
264  TFile *f = TFile::Open(filename, "RECREATE");
265 
266  savedir->cd();
267 
268  if (!f || f->IsZombie()) return 0;
269 
270  Event *event=new Event();
271  Event *ep = event;
272  TTree* eventtree= new TTree("EventTree", "Event Tree");
273  eventtree->SetDirectory(f);
274 
275  const Int_t buffersize=32000;
276  eventtree->Branch("event", "Event", &ep, buffersize, 1);
277  eventtree->AutoSave();
278 
279  Long64_t i=0;
280  Long64_t size_generated=0;
281 
282 // f->SetCompressionLevel(0); //no compression
283  Int_t ntrks = fNTracks;
284 
285  Info("GenerateFiles", "Generating %s", filename);
286  while (sizenevents--){
287  //event->Build(i++,fNTracksBench,0);
288  if (fNTracksMax > fNTracks) {
289  // Required to smear the number of tracks between [min,max]
290  ntrks = fNTracks + gRandom->Integer(fNTracksMax - fNTracks);
291  }
292  event->Build(i++, ntrks, 0);
293  size_generated+=eventtree->Fill();
294  }
295  nentries=eventtree->GetEntries();
296  Info("GenerateFiles", "%s generated with %lld entries", filename, nentries);
297  savedir = gDirectory;
298 
299  f = eventtree->GetCurrentFile();
300  f->cd();
301  eventtree->Write();
302  eventtree->SetDirectory(0);
303 
304  f->Close();
305  delete f;
306  f = 0;
307  eventtree->Delete();
308  event->Delete();
309  savedir->cd();
310 
311  return nentries;
312 }
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// The Process() function is called for each entry in the tree (or possibly
316 /// keyed object in the case of PROOF) to be processed. The entry argument
317 /// specifies which entry in the currently loaded tree is to be processed.
318 /// It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
319 /// to read either all or the required parts of the data. When processing
320 /// keyed objects with PROOF, the object is already loaded and is available
321 /// via the fObject pointer.
322 ///
323 /// This function should contain the "body" of the analysis. It can contain
324 /// simple or elaborate selection criteria, run algorithms on the data
325 /// of the event and typically fill histograms.
326 
327 Bool_t TSelEventGen::Process(Long64_t entry)
328 {
329  // WARNING when a selector is used with a TChain, you must use
330  // the pointer to the current TTree to call GetEntry(entry).
331  // The entry is always the local entry number in the current tree.
332  // Assuming that fChain is the pointer to the TChain being processed,
333  // use fChain->GetTree()->GetEntry(entry).
334 
335  TDSetElement *fCurrent = 0;
336  TPair *elemPair = 0;
337  if (fInput && (elemPair = dynamic_cast<TPair *>
338  (fInput->FindObject("PROOF_CurrentElement")))) {
339  if ((fCurrent = dynamic_cast<TDSetElement *>(elemPair->Value()))) {
340  Info("Process", "entry %lld: file: '%s'", entry, fCurrent->GetName());
341  } else {
342  Error("Process", "entry %lld: no file specified!", entry);
343  return kFALSE;
344  }
345  }
346 
347  // Generate
348  TString filename(fCurrent->GetName());
349  if (!fBaseDir.IsNull()) {
350  if (fBaseDir.Contains("<fn>")) {
351  filename = fBaseDir;
352  filename.ReplaceAll("<fn>", fCurrent->GetName());
353  } else {
354  filename.Form("%s/%s", fBaseDir.Data(), fCurrent->GetName());
355  }
356  }
357  TString fndset(filename);
358 
359  // Set the Url for remote access
360  TString seed = TString::Format("%s/%s", gSystem->HostName(), filename.Data()), dsrv;
361  TUrl basedirurl(filename, kTRUE);
362  if (!strcmp(basedirurl.GetProtocol(), "file")) {
363  TProofServ::GetLocalServer(dsrv);
364  TProofServ::FilterLocalroot(fndset, dsrv);
365  }
366 
367  //generate files
368  Long64_t neventstogenerate = fNEvents;
369 
370  Long64_t entries_file=0;
371  Long64_t filesize=0;
372  Bool_t filefound=kFALSE;
373  FileStat_t filestat;
374  TUUID uuid;
375  if (!fRegenerate && !gSystem->GetPathInfo(filename, filestat)) { //stat'ed
376  TFile *f = TFile::Open(filename);
377  if (f && !f->IsZombie()){
378  TTree* t = (TTree *) f->Get("EventTree");
379  if (t) {
380  entries_file = t->GetEntries();
381  if (entries_file == neventstogenerate) {
382  // File size seems to be correct, skip generation
383  Info("Process", "bench file (%s, entries=%lld) exists:"
384  " skipping generation.", filename.Data(), entries_file);
385  filesize = f->GetSize();
386  uuid = f->GetUUID();
387  filefound = kTRUE;
388  }
389  }
390  f->Close();
391  }
392  SafeDelete(f);
393  }
394 
395  // Make sure there is enough space left of the device, if local
396  TString bdir(fBaseDir);
397  bdir.ReplaceAll("<fn>", "");
398  if (!gSystem->AccessPathName(bdir)) {
399  Long_t devid, devbsz, devbtot, devbfree;
400  gSystem->GetFsInfo(bdir, &devid, &devbsz, &devbtot, &devbfree);
401  // Must be more than 10% of space and at least 1 GB
402  Long_t szneed = 1024 * 1024 * 1024, tomb = 1024 * 1024;
403  if (devbfree * devbsz < szneed || devbfree < 0.1 * devbtot) {
404  Error("Process", "not enough free space on device (%ld MB < {%ld, %ld} MB):"
405  " skipping generation of: %s",
406  (devbfree * devbsz) / tomb,
407  szneed / tomb, (Long_t) (0.1 * devbtot * devbsz / tomb),
408  filename.Data());
409  fStatus = TSelector::kAbortFile;
410  }
411  }
412 
413  if (!filefound) { // Generate
414  gRandom->SetSeed(static_cast<UInt_t>(TMath::Hash(seed)));
415  if (fGenerateFun) {
416  TString fargs = TString::Format("\"%s\",%lld", filename.Data(), neventstogenerate);
417  entries_file = (Long64_t) fGenerateFun->Exec(fargs);
418  } else {
419  entries_file = GenerateFiles(filename, neventstogenerate);
420  }
421 
422  TFile *f = TFile::Open(filename);
423  if (f && !f->IsZombie()) {
424  filesize = f->GetSize();
425  uuid = f->GetUUID();
426  f->Close();
427  } else {
428  Error("Process", "can not open generated file: %s", filename.Data());
429  fStatus = TSelector::kAbortFile;
430  return kFALSE;
431  }
432 
433  SafeDelete(f);
434  }
435 
436  // Add meta data to the file info
437  TFileInfoMeta* fimeta = new TFileInfoMeta("/EventTree", "TTree", entries_file);
438  TMD5* md5 = 0;
439  if (!strcmp(TUrl(filename,kTRUE).GetProtocol(), "file"))
440  md5 = TMD5::FileChecksum(filename);
441  TString md5s = (md5) ? md5->AsString() : "";
442  TFileInfo *fi = new TFileInfo(TString::Format("%s%s", dsrv.Data(), fndset.Data()),
443  filesize, uuid.AsString(), md5s.Data(), fimeta);
444  SafeDelete(md5);
445 
446  // Mark it as staged
447  fi->SetBit(TFileInfo::kStaged);
448 
449  // Add the fileinfo to the list
450  if (fFilesGenerated) fFilesGenerated->Add(fi);
451 
452  return kTRUE;
453 }
454 
455 ////////////////////////////////////////////////////////////////////////////////
456 /// The SlaveTerminate() function is called after all entries or objects
457 /// have been processed. When running with PROOF SlaveTerminate() is called
458 /// on each slave server
459 
460 void TSelEventGen::SlaveTerminate()
461 {
462  if (fFilesGenerated && fFilesGenerated->GetSize() > 0) {
463  fOutput->Add(fFilesGenerated);
464  Info("SlaveTerminate",
465  "list '%s' of files generated by this worker added to the output list",
466  fFilesGenerated->GetName());
467  } else {
468  if (!fFilesGenerated) {
469  Warning("SlaveTerminate", "no list of generated files defined!");
470  } else {
471  Warning("SlaveTerminate", "list of generated files is empty!");
472  }
473  }
474 }
475 
476 ////////////////////////////////////////////////////////////////////////////////
477 /// The Terminate() function is the last function to be called during
478 /// a query. It always runs on the client, it can be used to present
479 /// the results graphically or save the results to file.
480 
481 void TSelEventGen::Terminate()
482 {
483 }
484 
485 ////////////////////////////////////////////////////////////////////////////////
486 
487 void TSelEventGen::Print(Option_t *) const
488 {
489  Printf("fNEvents=%lld", fNEvents);
490  Printf("fBaseDir=%s", fBaseDir.Data());
491  Printf("fNTracks=%d", fNTracks);
492  Printf("fRegenerate=%d", fRegenerate);
493 }
494