StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
RunFastJet.C
1 //
2 // Pibero Djawotho <pibero@tamu.edu>
3 // Texas A&M University
4 // 31 Aug 2011
5 //
6 // $Id: RunFastJet.C,v 1.4 2012/07/16 21:50:26 pibero Exp $
7 //
8 // $Log: RunFastJet.C,v $
9 // Revision 1.4 2012/07/16 21:50:26 pibero
10 // Updated for easier reading...
11 //
12 // Revision 1.3 2012/03/10 23:10:36 pibero
13 // Added support for fastjet plugins
14 //
15 // Revision 1.2 2011/08/31 18:03:29 pibero
16 // Minor updates
17 //
18 //
19 
20 void RunFastJet(int nevents = 1e6,
21  const char* mudstfile = "root://xrdstar.rcf.bnl.gov:1095//home/starlib/home/starreco/reco/production2009_200Gev_Single/ReversedFullField/P10ic/2009/143/10143008/st_physics_10143008_raw_6020001.MuDst.root",
22  const char* jetfile = "jets.root",
23  const char* skimfile = "skim.root",
24  bool useL2 = false)
25 {
26  cout << "Read MuDst file:\t" << mudstfile << endl;
27  cout << "Write jet file:\t" << jetfile << endl;
28  cout << "Write skim file:\t" << skimfile << endl;
29 
30  gROOT->Macro("loadMuDst.C");
31  gROOT->Macro("LoadLogger.C");
32 
33  gSystem->Load("StDetectorDbMaker");
34  gSystem->Load("StTpcDb");
35  gSystem->Load("StDbUtilities");
36  gSystem->Load("StMcEvent");
37  gSystem->Load("StMcEventMaker");
38  gSystem->Load("StDaqLib");
39  gSystem->Load("StEmcRawMaker");
40  gSystem->Load("StEmcADCtoEMaker");
41  gSystem->Load("StEpcMaker");
42  gSystem->Load("StEmcSimulatorMaker");
43  gSystem->Load("StDbBroker");
44  gSystem->Load("St_db_Maker");
45  gSystem->Load("StEEmcUtil");
46  gSystem->Load("StEEmcDbMaker");
47  gSystem->Load("StSpinDbMaker");
48  gSystem->Load("StEmcTriggerMaker");
49  gSystem->Load("StTriggerUtilities");
50  gSystem->Load("StMCAsymMaker");
51  gSystem->Load("StRandomSelector");
52  gSystem->Load("libfastjet.so");
53  gSystem->Load("libCDFConesPlugin.so");
54  gSystem->Load("libEECambridgePlugin.so");
55  gSystem->Load("libJadePlugin.so");
56  gSystem->Load("libNestedDefsPlugin.so");
57  gSystem->Load("libSISConePlugin.so");
58  gSystem->Load("StJetFinder");
59  gSystem->Load("StJetSkimEvent");
60  gSystem->Load("StJets");
61  gSystem->Load("StJetEvent");
62  gSystem->Load("StJetMaker");
63  gSystem->Load("StTriggerFilterMaker");
64 
65  // The chain
66  StChain* chain = new StChain;
67 
68  // MuDst reader
69  StMuDstMaker* muDstMaker = new StMuDstMaker(0,0,"",mudstfile,"",100000,"MuDst");
70 
71  // MuDst DB
72  StMuDbReader* muDstDb = StMuDbReader::instance();
73 
74  // Trigger filter
75  StTriggerFilterMaker* filterMaker = new StTriggerFilterMaker;
76 
77  // 2009 pp500
78  filterMaker->addTrigger(230410); // JP1
79  filterMaker->addTrigger(230411); // JP2
80  filterMaker->addTrigger(230420); // AJP
81  filterMaker->addTrigger(230531); // BHT3
82 
83  // 2009 pp200
84  // http://www.star.bnl.gov/protected/common/common2009/trigger2009/triggers2009.html
85  // L2JetHigh
86  filterMaker->addTrigger(240650);
87  filterMaker->addTrigger(240651);
88  filterMaker->addTrigger(240652);
89  // JP1
90  filterMaker->addTrigger(240410);
91  filterMaker->addTrigger(240411);
92  // L2BGamma
93  filterMaker->addTrigger(240620);
94  // L2EGamma
95  filterMaker->addTrigger(240630);
96  filterMaker->addTrigger(240631);
97  // BHT3
98  filterMaker->addTrigger(240530);
99  // BBCMB-Cat2
100  filterMaker->addTrigger(240013);
101  filterMaker->addTrigger(240113);
102  filterMaker->addTrigger(240123);
103  filterMaker->addTrigger(240223);
104  // BBCMB-Cat3
105  filterMaker->addTrigger(240014);
106  filterMaker->addTrigger(240114);
107  filterMaker->addTrigger(240124);
108  filterMaker->addTrigger(240224);
109 
110  // star database
111  St_db_Maker* starDb = new St_db_Maker("StarDb","MySQL:StarDb");
112 
113  // Endcap database
114  StEEmcDbMaker* eemcDb = new StEEmcDbMaker;
115 
116  // Spin database
117  StSpinDbMaker* spinDb = new StSpinDbMaker;
118 
119  // Barrel ADC to energy maker
121  adc->saveAllStEvent(true);
122 
123  // Trigger simulator
124  StTriggerSimuMaker* simuTrig = new StTriggerSimuMaker;
125  simuTrig->useOnlineDB();
126  simuTrig->setMC(false); // Must be before individual detectors, to be passed
127  // BBC was not used in Run 9
128  //simuTrig->useBbc();
129  simuTrig->useBemc();
130  simuTrig->useEemc();
131  simuTrig->bemc->setConfig(StBemcTriggerSimu::kOffline);
132 
133  // L2 (only L2btowCalib, L2etowCalib, L2ped, L2jet in CVS as of 17 April 2010)
134  if (useL2) {
136  assert(simL2Mk);
137  simL2Mk->setSetupPath("/star/u/pibero/public/StarTrigSimuSetup/");
138  simL2Mk->setOutPath("./");
139  simuTrig->useL2(simL2Mk);
140  }
141 
142  // Skim event maker
143  StJetSkimEventMaker* skimEventMaker = new StJetSkimEventMaker("StJetSkimEventMaker",muDstMaker,skimfile);
144 
145  // Jet maker
146  StJetMaker2009* jetmaker = new StJetMaker2009;
147  jetmaker->setJetFile(jetfile);
148 
149  // Set analysis cuts for 12-point branch
150  StAnaPars* anapars12 = new StAnaPars;
151  anapars12->useTpc = true;
152  anapars12->useBemc = true;
153  anapars12->useEemc = true;
154 
155  // The classes available for correcting tower energy for tracks are:
156  // 1. StjTowerEnergyCorrectionForTracksMip
157  // 2. StjTowerEnergyCorrectionForTracksFraction
158  // 3. StjTowerEnergyCorrectionForTracksNull (default: no correction)
159  anapars12->setTowerEnergyCorrection(new StjTowerEnergyCorrectionForTracksFraction(1.00));
160 
161  // TPC cuts
162  anapars12->addTpcCut(new StjTrackCutFlag(0));
163  anapars12->addTpcCut(new StjTrackCutNHits(12));
164  anapars12->addTpcCut(new StjTrackCutPossibleHitRatio(0.51));
165  anapars12->addTpcCut(new StjTrackCutDca(3));
166  anapars12->addTpcCut(new StjTrackCutDcaPtDependent);
167  anapars12->addTpcCut(new StjTrackCutChi2(0,4));
168  anapars12->addTpcCut(new StjTrackCutPt(0.2,200));
169  anapars12->addTpcCut(new StjTrackCutEta(-2.5,2.5));
170  anapars12->addTpcCut(new StjTrackCutLastPoint(125));
171 
172  // BEMC cuts
173  anapars12->addBemcCut(new StjTowerEnergyCutBemcStatus(1));
174  anapars12->addBemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
175  anapars12->addBemcCut(new StjTowerEnergyCutEt(0.2));
176 
177  // EEMC cuts
178  anapars12->addEemcCut(new StjTowerEnergyCutBemcStatus(1));
179  anapars12->addEemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
180  anapars12->addEemcCut(new StjTowerEnergyCutEt(0.2));
181 
182  // Jet cuts
183  anapars12->addJetCut(new StProtoJetCutPt(5,200));
184  anapars12->addJetCut(new StProtoJetCutEta(-100,100));
185 
186  // Set analysis cuts for 5-point branch
187  StAnaPars* anapars5 = new StAnaPars;
188  anapars5->useTpc = true;
189  anapars5->useBemc = true;
190  anapars5->useEemc = true;
191 
192  // The classes available for correcting tower energy for tracks are:
193  // 1. StjTowerEnergyCorrectionForTracksMip
194  // 2. StjTowerEnergyCorrectionForTracksFraction
195  // 3. StjTowerEnergyCorrectionForTracksNull (default: no correction)
196  anapars5->setTowerEnergyCorrection(new StjTowerEnergyCorrectionForTracksFraction(1.00));
197 
198  // TPC cuts
199  anapars5->addTpcCut(new StjTrackCutFlag(0));
200  anapars5->addTpcCut(new StjTrackCutNHits(5));
201  anapars5->addTpcCut(new StjTrackCutPossibleHitRatio(0.51));
202  anapars5->addTpcCut(new StjTrackCutDca(3));
203  anapars5->addTpcCut(new StjTrackCutDcaPtDependent);
204  anapars5->addTpcCut(new StjTrackCutChi2(0,4));
205  anapars5->addTpcCut(new StjTrackCutPt(0.2,200));
206  anapars5->addTpcCut(new StjTrackCutEta(-2.5,2.5));
207 
208  // BEMC cuts
209  anapars5->addBemcCut(new StjTowerEnergyCutBemcStatus(1));
210  anapars5->addBemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
211  anapars5->addBemcCut(new StjTowerEnergyCutEt(0.2));
212 
213  // EEMC cuts
214  anapars5->addEemcCut(new StjTowerEnergyCutBemcStatus(1));
215  anapars5->addEemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
216  anapars5->addEemcCut(new StjTowerEnergyCutEt(0.2));
217 
218  // Jet cuts
219  anapars5->addJetCut(new StProtoJetCutPt(5,200));
220  anapars5->addJetCut(new StProtoJetCutEta(0.8,2.5));
221 
222  // Set analysis cuts for EMC branch
223  StAnaPars* anaparsEMC = new StAnaPars;
224  anaparsEMC->useTpc = true;
225  anaparsEMC->useBemc = true;
226  anaparsEMC->useEemc = true;
227 
228  // TPC cuts
229  anaparsEMC->addTpcCut(new StjTrackCutFlag(0));
230  anaparsEMC->addTpcCut(new StjTrackCutNHits(1000000));
231 
232  // BEMC cuts
233  anaparsEMC->addBemcCut(new StjTowerEnergyCutBemcStatus(1));
234  anaparsEMC->addBemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
235  anaparsEMC->addBemcCut(new StjTowerEnergyCutEt(0.2));
236 
237  // EEMC cuts
238  anaparsEMC->addEemcCut(new StjTowerEnergyCutBemcStatus(1));
239  anaparsEMC->addEemcCut(new StjTowerEnergyCutAdc(4,3)); // ADC-ped>4 AND ADC-ped>3*RMS
240  anaparsEMC->addEemcCut(new StjTowerEnergyCutEt(0.2));
241 
242  // Jet cuts
243  anaparsEMC->addJetCut(new StProtoJetCutPt(5,200));
244  anaparsEMC->addJetCut(new StProtoJetCutEta(-100,100));
245 
246  // Set anti-kt parameters
247  const double CONE_RADIUS = 0.6;
248 
249  StFastJetPars* AntiKtPars = new StFastJetPars;
250  AntiKtPars->setJetAlgorithm(StFastJetPars::antikt_algorithm);
251  AntiKtPars->setRparam(CONE_RADIUS);
252  AntiKtPars->setRecombinationScheme(StFastJetPars::E_scheme);
253  AntiKtPars->setStrategy(StFastJetPars::Best);
254  AntiKtPars->setPtMin(5.0);
255 
256  jetmaker->addBranch("ConeJets12",anapars12,AntiKtPars);
257  jetmaker->addBranch("ConeJets5",anapars5,AntiKtPars);
258  jetmaker->addBranch("ConeJetsEMC",anaparsEMC,AntiKtPars);
259 
260  chain->Init();
261  chain->EventLoop(nevents);
262 }
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.