StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
rdMuWana2011.C
1 /* How to use this macro.
2 To run w/o jets
3 
4 root4star -b -q 'rdMuWana2011.C(2e3,"st_W_12037041_raw_1400001.MuDst.root",0,0,5,7,"/star/data13/Magellan/reco/pp500_production_2011/ReversedFullField/P10k/2011/12037041/")'
5 
6 To produce jets
7 root4star -b -q 'rdMuWana2011.C(2e3,"st_W_12037041_raw_1400001.MuDst.root",0,1,5,7,"/star/data13/Magellan/reco/pp500_production_2011/ReversedFullField/P10k/2011/12037041/","/star/institutions/anl/balewski/jetTemp/")'
8 
9 */
10 
11 
12 
13 class StChain;
14 StChain *chain = 0;
15 
16 int spinSort = false;
17 bool isZ = false;
18 int geant = false;
19 
20 int rdMuWana2011(
21  int nEve = 1e6,
22  char* file = "st_W_12037063_raw_1380001_1201.MuDst.root",
23  int isMC = 0, // 0 = run9-data 200 = new MC w/ EEss in BFC
24  int useJetFinder = 0, // 0 - no jets = badWalgo; 1 generate jet trees; 2 read jet trees
25  int idL2BWtrg = 0, // offline Ids needed for real data
26  int idL2EWtrg = 0, // run 9 L2EW
27  // make those below empty for scheduler
28  char* muDir = "",
29  char* jetDir = "",
30  char* histDir = "",
31  char* wtreeDir = ""
32 )
33 {
34  char *eemcSetupPath = "/afs/rhic.bnl.gov/star/users/kocolosk/public/StarTrigSimuSetup/";
35 
36  if (isMC && useJetFinder == 2) geant = true;
37 
38  if (isMC) spinSort = false;
39 
40  string inputPathFile(file);
41 
42  size_t iLastSlash = inputPathFile.find_last_of("/");
43 
44  string inputPath("");
45 
46  if (iLastSlash != string::npos)
47  inputPath = inputPathFile.substr(0, iLastSlash);
48  else inputPath = "";
49 
50  string inputFile = inputPathFile.substr(iLastSlash + 1);
51 
52  printf("Input path: %s\n", inputPath.c_str());
53  printf("Input file: %s\n", inputFile.c_str());
54 
55  TString outF = inputFile;
56  outF = outF.ReplaceAll(".MuDst.root", "");
57  outF = outF.ReplaceAll(".lis", "");
58 
59  if (!isMC) {
60  //outF = file;
61  //outF = outF;
62  }
63  else { // new MC w/ working time stamp
64  //assert(2==5); M-C not unpacking not implemented
65  char *file1 = strstr(file, "cn100");
66  assert(file1);
67  file1--;
68  printf("file1: %s\n", file1);
69  outF = file1;
70  //outF.ReplaceAll(".MuDst.root","");
71  TString fileG = file;
72  fileG.ReplaceAll("MuDst", "geant");
73  }
74 
75  printf("Output file: %s\n", outF.Data());
76 
77  printf("TRIG ID: L2BW=%d, L2EW=%d isMC=%d useJetFinder=%d\n", idL2BWtrg, idL2EWtrg, isMC, useJetFinder );
78 
79  gROOT->LoadMacro("$STAR/StRoot/StMuDSTMaker/COMMON/macros/loadSharedLibraries.C");
80  loadSharedLibraries();
81  gROOT->LoadMacro("LoadLogger.C");
82  gMessMgr->SwitchOff("D");
83  gMessMgr->SwitchOn("I");
84 
85  assert( !gSystem->Load("StDaqLib"));
86 
87  // libraries below are needed for DB interface
88  assert( !gSystem->Load("StDetectorDbMaker")); // for St_tpcGasC
89  assert( !gSystem->Load("StTpcDb"));
90  assert( !gSystem->Load("StDbUtilities")); //for trigger simu
91  assert( !gSystem->Load("StDbBroker"));
92  assert( !gSystem->Load("St_db_Maker"));
93  assert( !gSystem->Load("StEEmcUtil"));
94  assert( !gSystem->Load("StEEmcDbMaker"));
95  assert( !gSystem->Load("StTriggerFilterMaker"));
96  assert( !gSystem->Load("StWalgo2011"));
97  assert( !gSystem->Load("StTriggerUtilities"));
98  assert( !gSystem->Load("StSpinDbMaker"));
99 
100  if (useJetFinder == 1 || useJetFinder == 2) { // jetfinder/jetreader libraries
101  cout << "BEGIN: loading jetfinder libs" << endl;
102  gSystem->Load("StEmcRawMaker");
103  gSystem->Load("StEmcADCtoEMaker");
104  gSystem->Load("StJetSkimEvent");
105  gSystem->Load("StJets");
106  gSystem->Load("StSpinDbMaker");
107  gSystem->Load("StEmcTriggerMaker");
108  gSystem->Load("StTriggerUtilities");
109  gSystem->Load("StMCAsymMaker");
110  gSystem->Load("StRandomSelector");
111  gSystem->Load("StJetEvent");
112  gSystem->Load("StJetFinder");
113  gSystem->Load("StJetMaker");
114  cout << "END: loading jetfinder libs" << endl;
115  }
116  else {
117  cout << "\nWARN: Jet are NOT read in, W-algo will not wrk properly\n " << endl;
118  }
119 
120  if (geant) {
121  // libraries for access to MC record
122  assert( !gSystem->Load("StMcEvent"));
123  assert( !gSystem->Load("StMcEventMaker"));
124 
125  // libraries for trigger simulator
126  assert( !gSystem->Load("StEmcSimulatorMaker"));
127  assert( !gSystem->Load("StEEmcSimulatorMaker"));
128  assert( !gSystem->Load("StEpcMaker"));
129  }
130 
131  cout << " loading done " << endl;
132 
133  // create chain
134  chain = new StChain("StChain");
135 
136  // create histogram storage array (everybody needs it):
137  TObjArray* HList = new TObjArray;
138 
139  if (geant) {
140  // get geant file
141  StIOMaker* ioMaker = new StIOMaker();
142  printf("\n %s \n\n", fileG.Data());
143  ioMaker->SetFile(fileG.Data());
144 
145  ioMaker->SetIOMode("r");
146  ioMaker->SetBranch("*", 0, "1"); //deactivate all branches
147  ioMaker->SetBranch("geantBranch", 0, "r"); //activate geant Branch
148  ioMaker->SetBranch("minimcBranch", 0, "r"); //activate geant Branch
149  }
150 
151  // Now we add Makers to the chain...
152  int maxFiles = 1000;
153 
154  StMuDstMaker *stMuDstMaker = new StMuDstMaker(0, 0, muDir, file, "MuDst.root", maxFiles);
155 
156  stMuDstMaker->SetStatus("*", 0);
157  stMuDstMaker->SetStatus("MuEvent", 1);
158  stMuDstMaker->SetStatus("EmcTow", 1);
159  stMuDstMaker->SetStatus("EmcSmde", 1);
160  stMuDstMaker->SetStatus("EmcSmdp", 1);
161  stMuDstMaker->SetStatus("PrimaryVertices", 1);
162  stMuDstMaker->SetStatus("GlobalTracks", 1);
163  stMuDstMaker->SetStatus("PrimaryTracks", 1);
164 
165  TChain* stMuDstMakerChain = stMuDstMaker->chain();
166 
167  assert(stMuDstMakerChain);
168 
169  int nEntries = (int) stMuDstMakerChain->GetEntries();
170 
171  if (nEntries < 0) {
172  Error("rdMuWana2011", "Invalid number of events %d", nEntries)
173  return -1;
174  }
175 
176  printf("Total number of events in muDst chain = %d\n", nEntries);
177 
178  //for EEMC, need full db access:
179  St_db_Maker *dbMk = new St_db_Maker("StarDb", "MySQL:StarDb", "MySQL:StarDb", "$STAR/StarDb");
180 
181  if (isMC == 0) {
182  // run 11 data
183  dbMk->SetFlavor("Wbose2", "bsmdeCalib"); // Willie's abs gains E-plane, run 9
184  dbMk->SetFlavor("Wbose2", "bsmdpCalib"); // P-plane
185  dbMk->SetFlavor("sim", "bemcCalib"); // use ideal gains for real data
186  dbMk->SetFlavor("sim", "eemcPMTcal"); // use ideal gains for 2011 real data as well
187  }
188  else {
189  printf("???? unforeseen MC flag, ABORT\n");
190  assert(1 == 2);
191  }
192 
193 
194  //.... load EEMC database
195  StEEmcDbMaker* mEEmcDatabase = new StEEmcDbMaker("eemcDb");
196 
197 #if 0 // drop abs lumi for now
198  if (!isMC && strstr(file, "fillListPhys")) {
199  StTriggerFilterMaker *filterMaker = new StTriggerFilterMaker;
200  filterMaker->addTrigger(230420); // AJP
201  filterMaker->addTrigger(230411); // JP2
202  filterMaker->addTrigger(bht3ID); // regular W -> e+ analysis
203  }
204 #endif
205 
206  if (geant) {
207  StMcEventMaker *mcEventMaker = new StMcEventMaker();
208  mcEventMaker->doPrintEventInfo = false;
209  mcEventMaker->doPrintMemoryInfo = false;
210  }
211 
212  if (geant && useJetFinder != 1) { // only use trigger simulator in W algo
213  //BEMC simulator:
214  StEmcSimulatorMaker* emcSim = new StEmcSimulatorMaker(); //use this instead to "redo" converstion from geant->adc
215  emcSim->setCalibSpread(kBarrelEmcTowerId, 0.15); //spread gains by 15%
216  StEmcADCtoEMaker *bemcAdc = new StEmcADCtoEMaker();//for real data this sets calibration and status
217 
218  //EEMC simulator:
219  new StEEmcDbMaker("eemcDb");
220  StEEmcSlowMaker *slowSim = new StEEmcSlowMaker("slowSim");
221  //slowSim->setSamplingFraction(0.0384); // effectively scales all Tower energies with a factor of 1.3 (for old private filtered simu only!)
222  slowSim->setAddPed(true);
223  slowSim->setSmearPed(true);
224 
225  //Get TriggerMaker
226  StTriggerSimuMaker *simuTrig = new StTriggerSimuMaker("StarTrigSimu");
227  simuTrig->setHList(HList);
228  simuTrig->setMC(isMC); // must be before individual detectors, to be passed
229  simuTrig->useBbc();
230  simuTrig->useEemc(0);//default=0:just process ADC, 1,2:comp w/trgData,see .
231  simuTrig->eemc->setSetupPath(eemcSetupPath);
232  simuTrig->useBemc();
233  simuTrig->bemc->setConfig(2);
234  }
235 
236  TString jetFile = jetDir;
237 
238  //.... Jet finder code ....
239  if (useJetFinder > 0) {
240  //ds TString jetFile = jetDir;
241  //ds jetFile+="jets_"+outF+".root";
242  jetFile += "jets_";
243  jetFile += outF + ".root";
244  cout << "BEGIN: running jet finder/reader =" << jetFile << "=" << endl;
245  }
246 
247 
248  if (useJetFinder == 1)
249  { //{{{ // run jet finder
250  double pi = atan(1.0) * 4.0;
251  // Makers for clusterfinding
252  StSpinDbMaker *uspDbMaker = new StSpinDbMaker("spinDb");
253  StEmcADCtoEMaker *adc = new StEmcADCtoEMaker();
254 
255  //here we also tag whether or not to do the swap:
256  bool doTowerSwapFix = true;
257  bool use2003TowerCuts = false;
258  bool use2006TowerCuts = true;
259 
260  //4p maker using 100% tower energy correction
261  StBET4pMaker* bet4pMakerFrac100 = new StBET4pMaker("BET4pMakerFrac100", stMuDstMaker, doTowerSwapFix, new StjTowerEnergyCorrectionForTracksFraction(1.0));
262  bet4pMakerFrac100->setUse2003Cuts(use2003TowerCuts);
263  bet4pMakerFrac100->setUseEndcap(true);
264  bet4pMakerFrac100->setUse2006Cuts(use2006TowerCuts);
265 
266  //4p maker using 100% tower energy correction (no endcap)
267  StBET4pMaker* bet4pMakerFrac100_noEEMC = new StBET4pMaker("BET4pMakerFrac100_noEEMC", stMuDstMaker, doTowerSwapFix, new StjTowerEnergyCorrectionForTracksFraction(1.0));
268  bet4pMakerFrac100_noEEMC->setUse2003Cuts(use2003TowerCuts);
269  bet4pMakerFrac100_noEEMC->setUseEndcap(false);
270  bet4pMakerFrac100_noEEMC->setUse2006Cuts(use2006TowerCuts);
271 
272  //Instantiate the JetMaker and SkimEventMaker
273  StJetMaker* emcJetMaker = new StJetMaker("emcJetMaker", stMuDstMaker, jetFile);
274  //StJetSkimEventMaker* skimEventMaker = new StJetSkimEventMaker("StJetSkimEventMaker", stMuDstMaker,outSkimFile);
275 
276  //set the analysis cuts: (see StJetMaker/StppJetAnalyzer.h -> class StppAnaPars )
277  StppAnaPars* anapars = new StppAnaPars();
278  anapars->setFlagMin(0); //track->flag() > 0
279  anapars->setNhits(12); //track->nHitsFit()>12
280  anapars->setCutPtMin(0.2); //track->pt() > 0.2
281  anapars->setAbsEtaMax(2.0); //abs(track->eta())<1.6
282  anapars->setJetPtMin(3.5);
283  anapars->setJetEtaMax(100.0);
284  anapars->setJetEtaMin(0);
285  anapars->setJetNmin(0);
286 
287  //Setup the cone finder (See StJetFinder/StConeJetFinder.h -> class StConePars)
288  StConePars* cpars = new StConePars();
289  cpars->setGridSpacing(105, -3.0, 3.0, 120, -pi, pi); //include EEMC
290  cpars->setConeRadius(0.7); // default=0.7
291  cpars->setSeedEtMin(0.5);
292  cpars->setAssocEtMin(0.1);
293  cpars->setSplitFraction(0.5); //default=0.5. if 0.3 less split?
294  cpars->setPerformMinimization(true);
295  cpars->setAddMidpoints(true);
296  cpars->setRequireStableMidpoints(true);
297  cpars->setDoSplitMerge(true);
298 
299  cpars->setDebug(false);
300 
301  emcJetMaker->addAnalyzer(anapars, cpars, bet4pMakerFrac100, "ConeJets12_100"); //100% subtraction
302  emcJetMaker->addAnalyzer(anapars, cpars, bet4pMakerFrac100_noEEMC, "ConeJets12_100_noEEMC"); //100% subtraction (no Endcap)
303 
304  chain->Init();
305  chain->ls(3);
306 
307  char txt[1000];
308  //---------------------------------------------------
309  int eventCounter = 0;
310  int t1 = time(0);
311  TStopwatch tt;
312 
313  for (Int_t iev = 0; iev < nEntries; iev++) {
314  if (eventCounter >= nEve) break;
315  chain->Clear();
316  int stat = chain->Make();
317  if (stat != kStOk && stat != kStSkip) break; // EOF or input error
318  eventCounter++;
319  }
320 
321  cout << "run " << file << " nEve=" << eventCounter << " total ";
322  tt.Print();
323  printf("****************************************** \n");
324 
325  int t2 = time(0);
326  if (t2 == t1) t2 = t1 + 1;
327  float tMnt = (t2 - t1) / 60.;
328  float rate = 1.*eventCounter / (t2 - t1);
329 
330  printf("jets sorting done %d of nEve= %d, CPU rate= %.1f Hz, total time %.1f minute(s) \n\n", eventCounter, nEntries, rate, tMnt);
331 
332  cout << "END: jet finder " << endl;
333 
334  return 1;
335  } //}}}
336 
337  if (useJetFinder == 2) {
338  cout << "Configure to read jet trees " << endl;
339  StJetReader *jetReader = new StJetReader;
340  jetReader->InitFile(jetFile);
341  }
342 
343  //.... W reconstruction code ....
344  St2011WMaker *WmuMk = new St2011WMaker();
345 
346  if (isMC) { // MC specific
347  WmuMk->setMC(isMC); //pass "version" of MC to maker
348  //WmuMk->setJetNeutScaleMC(1.0);
349  //WmuMk->setJetChrgScaleMC(1.0);
350  }
351  else { // real data
352  WmuMk->setTrigID(idL2BWtrg, idL2EWtrg);
353  }
354 
355  TString treeFileName = wtreeDir;
356  treeFileName += outF;
357  treeFileName += ".Wtree.root";
358 
359  WmuMk->setTreeName(treeFileName);
360 
361  if (useJetFinder == 2)
362  WmuMk->setJetTreeBranch("ConeJets12_100", "ConeJets12_100_noEEMC"); //select jet tree braches used
363 
364  WmuMk->setMaxDisplayEve(100); // only first N events will get displayed
365  //set energy scale (works for data and MC - be careful!)
366  //WmuMk->setBtowScale(1.0);
367  //WmuMk->setEtowScale(1.0);
368 
369  /* evaluation of result, has full acess to W-algo internal data
370  including overwrite - be careful */
371 
372  St2011pubWanaMaker *WpubMk = new St2011pubWanaMaker("pubJan");
373  WpubMk->attachWalgoMaker(WmuMk);
374 
375  //Collect all output histograms
376  //already defined this above: TObjArray* HList=new TObjArray;
377  WmuMk->setHList(HList);
378  WpubMk->setHList(HList);
379 
380  StSpinDbMaker *spDb = 0;
381  if (spinSort) {
382  spDb = new StSpinDbMaker("spinDb");
383  enum {mxSM = 5}; // to study eta-cuts, drop Q/PT cut
384  St2011pubSpinMaker *spinMkA[mxSM];
385 
386  for (int kk = 0; kk < mxSM; kk++) {
387  char ttx[100];
388  sprintf(ttx, "%cspin", 'A' + kk);
389  printf("add spinMaker %s %d \n", ttx, kk);
390  spinMkA[kk] = new St2011pubSpinMaker(ttx);
391  spinMkA[kk]->attachWalgoMaker(WmuMk);
392  spinMkA[kk]->attachSpinDb(spDb);
393  spinMkA[kk]->setHList(HList);
394  if (kk == 1) spinMkA[kk]->setEta(-1., 0.);
395  if (kk == 2) spinMkA[kk]->setEta(0, 1.);
396  if (kk == 3) spinMkA[kk]->setQPT(-1); // disable Q/PT cut
397  if (kk == 4) spinMkA[kk]->setNoEEMC();
398  }
399  }
400 
401  if (geant) {
402  pubMcMk = new St2009pubMcMaker("pubMc");
403  pubMcMk->attachWalgoMaker(WmuMk);
404  pubMcMk->setHList(HList);
405  }
406 
407  if (isZ) {
408  ZMk = new St2011ZMaker("Z");
409  ZMk->attachWalgoMaker(WmuMk);
410  ZMk->setHList(HList);
411  ZMk->setNearEtFrac(0.88);
412  ZMk->setClusterMinEt(15);
413  ZMk->setPhi12Min(3.1416 / 2.);
414  ZMk->setMinZMass(73.); // Zmass -20%
415  ZMk->setMaxZMass(114.);// Zmass +20%
416  }
417 
418  chain->Init();
419  chain->ls(3);
420 
421  char txt[1000];
422  int eventCounter = 0;
423  int t1 = time(0);
424 
425  TStopwatch tt;
426 
427  for (Int_t iev = 0; iev < nEntries; iev++) {
428  Info("rdMuWana2011", "Analyzing event %d", iev);
429  if (eventCounter >= nEve) break;
430  chain->Clear();
431  int stat = chain->Make();
432  if (stat != kStOk && stat != kStSkip) break; // EOF or input error
433  eventCounter++;
434  }
435 
436  //chain->Finish();
437 
438  cout << "run " << file << " nEve=" << eventCounter << " total ";
439  tt.Print();
440 
441  printf("****************************************** \n");
442 
443  int t2 = time(0);
444  if (t2 == t1) t2 = t1 + 1;
445  float tMnt = (t2 - t1) / 60.;
446  float rate = 1.*eventCounter / (t2 - t1);
447 
448  printf("#sorting %s done %d of nEve= %d, CPU rate= %.1f Hz, total time %.1f minute(s) \n\n", file, eventCounter, nEntries, rate, tMnt);
449 
450 
451  TString histFileName = histDir;
452  histFileName += outF;
453  histFileName += ".wana.hist.root";
454 
455  cout << "Output histo file " << histFileName << endl;
456 
457  hf = new TFile(histFileName, "recreate");
458 
459  if (hf->IsOpen()) {
460  //HList->ls();
461  HList->Write();
462  printf("\n Histo saved -->%s<\n", outFh.Data());
463  }
464  else {
465  printf("\n Failed to open Histo-file -->%s<, continue\n", outFh.Data());
466  }
467 
468  //WmuMk->Finish();
469 
470  return 2;
471 }
472 
473 
474 // $Log: rdMuWana2011.C,v $
475 // Revision 1.5 2012/08/21 17:40:15 stevens4
476 // Revert to previous version
477 //
478 // Revision 1.3 2012/03/19 23:45:23 smirnovd
479 // Major clean up. Removed hardcoded references etc.
480 //
481 // Revision 1.2 2012/03/12 23:11:30 smirnovd
482 // *** empty log message ***
483 //
484 // Revision 1.1 2011/02/10 20:33:35 balewski
485 // start
486 //
void setConeRadius(double v)
Set cone radius:
Definition: StConePars.h:65
Definition: FJcore.h:367
gathers all results from W-analysis, Jan&#39;s analysis
maker to retrieve info from geant.root files for comparison with reco quantities from MC ...
void setAddMidpoints(bool v)
Add seeds at midpoints?
Definition: StConePars.h:56
virtual void SetIOMode(Option_t *iomode="w")
number of transactions
Definition: StIOInterFace.h:35
virtual void InitFile(const char *file)
Recover the TTree from file and prepare for reading.
Definition: StJetReader.cxx:61
void setAddPed(Bool_t a=true)
Add pedestal offsets from DB.
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
void setSplitFraction(double v)
split jets if E_shared/E_neighbor&gt;splitFraction
Definition: StConePars.h:50
void setGridSpacing(int nEta, double etaMin, double etaMax, int nPhi, double phiMin, double phiMax)
Set the grid spacing:
Definition: StConePars.h:36
void setDoSplitMerge(bool v)
Do Split/Merge step?
Definition: StConePars.h:59
spin sorting of W&#39;s
void setDebug(bool v)
Toggle debug streams on/off.
Definition: StConePars.h:68
Filling of all StMcEvent classes from g2t tables Transform all the data in the g2t tables into the co...
void setSeedEtMin(double v)
minimum et threshold to be considered a seed
Definition: StConePars.h:44
Definition: Stypes.h:49
virtual void ls(Option_t *option="") const
Definition: TDataSet.cxx:495
Bool_t doPrintMemoryInfo
lots of screen output
virtual Int_t Make()
Definition: StChain.cxx:110
void setAssocEtMin(double v)
minimum et threshold to be considered for addition to the seed
Definition: StConePars.h:47
TChain * chain()
In read mode, returns pointer to the chain of .MuDst.root files that where selected.
Definition: StMuDstMaker.h:426
void SetStatus(const char *arrType, int status)
void setSmearPed(Bool_t s=true)
Smear the pedestal with sigma from DB.
muDst based extraction of W-signal from pp500 data from 2011
Definition: St2011WMaker.h:49
void setRequireStableMidpoints(bool v)
Require stable midpoints?
Definition: StConePars.h:62
Slow simulator for EEMC.
Definition: Stypes.h:41
uses tree from W-algo to find Zs
Definition: St2011ZMaker.h:26
void setCalibSpread(StDetectorId det, float spread)
smear simulator calibration coefficients using Gaussian with this RMS.
void setPerformMinimization(bool v)
Let jet wander to minimum?
Definition: StConePars.h:53