StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
RunJetFinder2009pro_ue.C
1 //
2 // Grant Webb <gdwebb@bnl.gov>
3 // Brookhaven National Laboratory
4 // 10 August 2015
5 // Added UE to the JetFinder
6 //
7 
8 void RunJetFinder2009pro_ue(int nevents = 1E3,
9  const char* mudstfile = "/star/institutions/uky/gdwebb/UE_histos/st_physics_10103041_raw_8030001.MuDst.root",
10  const char* jetfile = "jets_copy.root",
11  const char* skimfile = "skim_copy.root",
12  const char* uefile = "ueTree_copy.root",
13  int mEmbed = 0,
14  int mPythia = 0,
15  bool useL2 = false)
16 {
17  cout << "nevents = " << nevents << endl;
18  cout << "mudstfile = " << mudstfile << endl;
19  cout << "jetfile = " << jetfile << endl;
20  cout << "skimfile = " << skimfile << endl;
21  cout << "uefile = " << uefile << endl;
22 
23  gROOT->Macro("loadMuDst.C");
24  gROOT->Macro("LoadLogger.C");
25 
26  gSystem->Load("StDetectorDbMaker");
27  gSystem->Load("StTpcDb");
28  gSystem->Load("StDbUtilities");
29  gSystem->Load("StMcEvent");
30  gSystem->Load("StMcEventMaker");
31  gSystem->Load("StDaqLib");
32  gSystem->Load("StEmcRawMaker");
33  gSystem->Load("StEmcADCtoEMaker");
34  gSystem->Load("StEpcMaker");
35  gSystem->Load("StEmcSimulatorMaker");
36  gSystem->Load("StDbBroker");
37  gSystem->Load("St_db_Maker");
38  gSystem->Load("StEEmcUtil");
39  gSystem->Load("StEEmcDbMaker");
40  gSystem->Load("StSpinDbMaker");
41  gSystem->Load("StEmcTriggerMaker");
42  gSystem->Load("StTriggerUtilities");
43  gSystem->Load("StMCAsymMaker");
44  gSystem->Load("StRandomSelector");
45 
46  gSystem->Load("libfastjet.so");
47  gSystem->Load("libsiscone.so");
48  gSystem->Load("libsiscone_spherical.so");
49  gSystem->Load("libfastjetplugins.so");
50 
51  // gSystem->Load("libfastjet.so");
52  // gSystem->Load("libCDFConesPlugin.so");
53  // gSystem->Load("libEECambridgePlugin.so");
54  // gSystem->Load("libJadePlugin.so");
55  // gSystem->Load("libNestedDefsPlugin.so");
56  // gSystem->Load("libSISConePlugin.so");
57  gSystem->Load("StJetFinder");
58  gSystem->Load("StJetSkimEvent");
59  gSystem->Load("StJets");
60  gSystem->Load("StJetEvent");
61  gSystem->Load("StUeEvent");
62  gSystem->Load("StJetMaker");
63  gSystem->Load("StTriggerFilterMaker");
64 
65  StChain* chain = new StChain;
66 
67  // MuDst reader
68  StMuDstMaker* muDstMaker = new StMuDstMaker(0,0,"",mudstfile,"",100000,"MuDst");
69 
70  // MuDst DB
71  StMuDbReader* muDstDb = StMuDbReader::instance();
72 
73  // Trigger filter
74  StTriggerFilterMaker* filterMaker = new StTriggerFilterMaker;
75 
76  // 2009 pp500
77  filterMaker->addTrigger(230410); // JP1
78  filterMaker->addTrigger(230411); // JP2
79  filterMaker->addTrigger(230420); // AJP
80  filterMaker->addTrigger(230531); // BHT3
81 
82  // star database
83  St_db_Maker* starDb = new St_db_Maker("StarDb","MySQL:StarDb");
84 
85  // Endcap database
86  StEEmcDbMaker* eemcDb = new StEEmcDbMaker;
87 
88  // Spin database
89  StSpinDbMaker* spinDb = new StSpinDbMaker;
90 
91  // Barrel ADC to energy maker
93  adc->saveAllStEvent(true);
94 
95  // Trigger simulator
96  StTriggerSimuMaker* simuTrig = new StTriggerSimuMaker;
97  simuTrig->useOnlineDB(); // for trigger definitions and thresholds
98  simuTrig->setMC(false); // Must be before individual detectors, to be passed
99 
100  simuTrig->useBemc();
101  simuTrig->useEemc();
102  simuTrig->bemc->setConfig(StBemcTriggerSimu::kOffline);
103 
104  // L2 (only L2btowCalib, L2etowCalib, L2ped, L2jet in CVS as of 17 April 2010)
105  if (useL2) {
107  assert(simL2Mk);
108  simL2Mk->setSetupPath("/star/u/pibero/public/StarTrigSimuSetup/");
109  simL2Mk->setOutPath("./");
110  simuTrig->useL2(simL2Mk);
111  }
112 
113  // Skim event maker
114  StJetSkimEventMaker* skimEventMaker = new StJetSkimEventMaker("StJetSkimEventMaker",muDstMaker,skimfile);
115 
116  // Jet maker
117  StJetMaker2009* jetmaker = new StJetMaker2009;
118  jetmaker->setJetFile(jetfile);
119  // UE maker
120  StUEMaker2009* uemaker = new StUEMaker2009;
121  uemaker->setUeFile(uefile);
122 
123  // Set analysis cuts for 12-point branch
124  StAnaPars* anapars12 = new StAnaPars;
125  anapars12->useTpc = true;
126  anapars12->useBemc = true;
127  anapars12->useEemc = false;
128 
129  // The classes available for correcting tower energy for tracks are:
130  // 1. StjTowerEnergyCorrectionForTracksMip
131  // 2. StjTowerEnergyCorrectionForTracksFraction
132  // 3. StjTowerEnergyCorrectionForTracksNull (default: no correction)
133  anapars12->setTowerEnergyCorrection(new StjTowerEnergyCorrectionForTracksFraction(1.00));
134 
135  //
136  // 1. StjTrackPtFraction
137  // 2. StjTowerEnergyFraction
138  // The input parameter is the fraction which you want to increase or decrease the track/tower pT/energy.
139  // i.e. If you want to decrease the track pT by 5%, you should put in -0.05 as the input parameter for StjTrackPtFraction
140  // anapars12->setTrackShift(new StjTrackPtFraction(-0.05));
141  // anapars12->setTowerShift(new StjTowerEnergyFraction(-0.05));
142 
143  // TPC cuts
144  anapars12->addTpcCut(new StjTrackCutFlag(0));
145  anapars12->addTpcCut(new StjTrackCutNHits(12));
146  anapars12->addTpcCut(new StjTrackCutPossibleHitRatio(0.51));
147  anapars12->addTpcCut(new StjTrackCutDca(3));
148  anapars12->addTpcCut(new StjTrackCutTdcaPtDependent);
149  anapars12->addTpcCut(new StjTrackCutPt(0.2,200));
150  anapars12->addTpcCut(new StjTrackCutEta(-2.5,2.5));
151  anapars12->addTpcCut(new StjTrackCutLastPoint(125));
152 
153  // BEMC cuts
154  anapars12->addBemcCut(new StjTowerEnergyCutBemcStatus(1));
155  anapars12->addBemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
156  anapars12->addBemcCut(new StjTowerEnergyCutEt(0.2));
157 
158  // EEMC cuts
159  anapars12->addEemcCut(new StjTowerEnergyCutBemcStatus(1));
160  anapars12->addEemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
161  anapars12->addEemcCut(new StjTowerEnergyCutEt(0.2));
162 
163  // Jet cuts
164  anapars12->addJetCut(new StProtoJetCutPt(5,200));
165  anapars12->addJetCut(new StProtoJetCutEta(-100,100));
166 
167  //---------Region Criteria---------
168  anapars12_toward = anapars12; // Toward Region for Tracks and Towers
169  anapars12_toward->setTrackRegion(new StjTrackRegion(60.0,-60.0,1.0));
170  anapars12_toward->setTowerRegion(new StjTowerRegion(60.0,-60.0,1.0));
171 
172  anapars12_away = anapars12; // Away Region for Tracks and Towers
173  anapars12_away->setTrackRegion(new StjTrackRegion(120.0,-120.0,1.0));
174  anapars12_away->setTowerRegion(new StjTowerRegion(120.0,-120.0,1.0));
175 
176  anapars12_transPlus = anapars12; // Trans Plus for Tracks and Towers
177  anapars12_transPlus->setTrackRegion(new StjTrackRegion(120.0,60.0,1.0));
178  anapars12_transPlus->setTowerRegion(new StjTowerRegion(120.0,60.0,1.0));
179 
180  anapars12_transMinus = anapars12; // Trans Minus for Tracks and Towers
181  anapars12_transMinus->setTrackRegion(new StjTrackRegion(-60.0,-120.0,1.0));
182  anapars12_transMinus->setTowerRegion(new StjTowerRegion(-60.0,-120.0,1.0));
183  //---------------------------------
184 
185  //Jet Area
186  StFastJetAreaPars *JetAreaPars = new StFastJetAreaPars;
187  JetAreaPars->setGhostArea(0.04);
188  // Set anti-kt R=0.6 parameters
189  StFastJetPars* AntiKtR060Pars = new StFastJetPars;
190  AntiKtR060Pars->setJetAlgorithm(StFastJetPars::antikt_algorithm);
191  AntiKtR060Pars->setRparam(0.6);
192  AntiKtR060Pars->setRecombinationScheme(StFastJetPars::E_scheme);
193  AntiKtR060Pars->setStrategy(StFastJetPars::Best);
194  AntiKtR060Pars->setPtMin(5.0);
195  AntiKtR060Pars->setJetArea(JetAreaPars);
196 
197  // Set anti-kt R=0.4 parameters
198  StFastJetPars* AntiKtR040Pars = new StFastJetPars;
199  AntiKtR040Pars->setJetAlgorithm(StFastJetPars::antikt_algorithm);
200  AntiKtR040Pars->setRparam(0.4);
201  AntiKtR040Pars->setRecombinationScheme(StFastJetPars::E_scheme);
202  AntiKtR040Pars->setStrategy(StFastJetPars::Best);
203  AntiKtR040Pars->setPtMin(5.0);
204  AntiKtR040Pars->setJetArea(JetAreaPars);
205 
206  // Set anti-kt R=0.5 parameters
207  StFastJetPars* AntiKtR050Pars = new StFastJetPars;
208  AntiKtR050Pars->setJetAlgorithm(StFastJetPars::antikt_algorithm);
209  AntiKtR050Pars->setRparam(0.5);
210  AntiKtR050Pars->setRecombinationScheme(StFastJetPars::E_scheme);
211  AntiKtR050Pars->setStrategy(StFastJetPars::Best);
212  AntiKtR050Pars->setPtMin(5.0);
213  AntiKtR050Pars->setJetArea(JetAreaPars);
214 
215  jetmaker->addBranch("AntiKtR060NHits12",anapars12,AntiKtR060Pars);
216  jetmaker->addBranch("AntiKtR040NHits12",anapars12,AntiKtR040Pars);
217  jetmaker->addBranch("AntiKtR050NHits12",anapars12,AntiKtR050Pars);
218 
219  uemaker->addBranch("toward",anapars12_toward,"AntiKtR060NHits12");
220  uemaker->addBranch("away",anapars12_away,"AntiKtR060NHits12");
221  uemaker->addBranch("transP",anapars12_transPlus,"AntiKtR060NHits12");
222  uemaker->addBranch("transM",anapars12_transMinus,"AntiKtR060NHits12");
223 
224 
225  // Run
226  chain->Init();
227  chain->EventLoop(nevents);
228 }
229 
static const int Best
automatic selection of the best (based on N)
void useOnlineDB()
Choose DB to access trigger definitions and thresholds.
static const int E_scheme
Definition: StFastJetPars.h:89
static const int antikt_algorithm
Definition: StFastJetPars.h:65
void saveAllStEvent(Bool_t a)
Set to kTRUE if all hits are to be saved on StEvent.