StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
St2009pubSpinMaker.cxx
1 // $Id: St2009pubSpinMaker.cxx,v 1.12 2010/05/01 01:31:45 balewski Exp $
2 //
3 //*-- Author : Jan Balewski, MIT
4 //
5 #include <StMuDSTMaker/COMMON/StMuDstMaker.h>
6 #include <StMuDSTMaker/COMMON/StMuDst.h>
7 #include <StMuDSTMaker/COMMON/StMuEvent.h>
8 #include <StEvent/StEventInfo.h> // just to get time
9 #include "StSpinPool/StSpinDbMaker/StSpinDbMaker.h"
10 
11 #include "St2009WMaker.h"
12 
13 #include "St2009pubSpinMaker.h"
14 
15 ClassImp(St2009pubSpinMaker)
16 
17 //_____________________________________________________________________________
18 //
19 St2009pubSpinMaker::St2009pubSpinMaker(const char *name):StMaker(name){
20  wMK=0;HList=0;
21  core=name;
22  par_QPTlow=0.010;
23 
24  par_QPThighET0=25;
25  par_QPThighET1=50;
26  par_QPThighA=0.08;
27  par_QPThighB=0.0013;
28  par_leptonEta1=-1.; par_leptonEta2=1.;
29  par_useNoEEMC=0;
30  }
31 
32 
33 //_____________________________________________________________________________
34 //
35 Int_t
36 St2009pubSpinMaker::Init(){
37  assert(wMK);
38  assert(HList);
39  initHistos();
40  nRun=0;
41  hbxIdeal=0;
42  return StMaker::Init();
43 }
44 
45 
46 //_____________________________________________________________________________
47 //
48 Int_t
49 St2009pubSpinMaker::FinishRun (int runNo){
50  char txt[1000];
51  sprintf(txt,"events T= %d %d",Tfirst,Tlast);
52  printf("Finish run=%d , events time range %s\n",runNo,txt);
53  hbxIdeal->GetYaxis()->SetTitle(txt);
54  return kStOK;
55 }
56 
57 //_____________________________________________________________________________
58 //
59 Int_t
60 St2009pubSpinMaker::InitRun (int runNo){
61  assert(runNo>= 10081007); // F10407, first pp500 long fill with defined pol pattern
62  assert(runNo<=10103046); // F10536, last pp500 fill in run9
63 
64  char txt[1000],txt0[100];
65  sprintf(txt0,"bxIdeal%d",nRun);
66  sprintf(txt,"intended fill pattern R%d-%d vs. bXing; %s", runNo,nRun,spinDb->getV124comment());
67  nRun++;
68  Tfirst=int(2e9); Tlast=-Tfirst;
69  hbxIdeal=new TH1F(core+txt0,txt,128,-0.5,127.5);
70  hbxIdeal->SetFillColor(kYellow);
71  HList->Add(hbxIdeal);
72 
73  spinDb->print(0); // 0=short, 1=huge
74  for(int bx=0;bx<120;bx++){
75  if(spinDb->isBXfilledUsingInternalBX(bx)) hbxIdeal->Fill(bx);
76  }
77 
78  sprintf(txt,"bXing= bx48+off=%d",spinDb->BX48offset());
79  hA[3]->GetXaxis()->SetTitle(txt);
80 
81  sprintf(txt,"bXing= bx7+off=%d",spinDb->BX7offset());
82  hA[4]->GetXaxis()->SetTitle(txt);
83 
84  LOG_INFO<<Form("::InitRun(%d) done, W-spin sorting params: exclude |Q/PT| < %.2f OR |Q/PT| above line %.3f*(ET-%.1f)-%6e if ET<%.1f, for AL use leptonEta in[%.1f,%.1f] useNoEEMC=%d", runNo,
85  par_QPTlow,par_QPThighET0, par_QPThighA ,par_QPThighB,par_QPThighET1,par_leptonEta1, par_leptonEta2,par_useNoEEMC
86  )<<endm;
87  return kStOK;
88 }
89 
90 //_____________________________________________________________________________
91 //
92 Int_t
94  int T=wMK->mMuDstMaker->muDst()->event()->eventInfo().time();
95  if(Tlast<T) Tlast=T;
96  if(Tfirst>T) Tfirst=T;
97 
98  bXingSort();
99  return kStOK;
100 }
101 
102 //_____________________________________________________________________________
103 //
104 void
105 St2009pubSpinMaker::bXingSort(){
106  //has access to whole W-algo-maker data via pointer 'wMK'
107 
108  hA[0]->Fill("inp",1.);
109 
110  assert(spinDb->isValid()); // all 3 DB records exist
111  assert(spinDb->isPolDirLong()); // you do not want mix Long & Trans by accident
112 
113  if(wMK->wEve.vertex.size()<=0) return;
114  //......... require: L2W-trig (ET or rnd) & vertex is reasonable .......
115 
116  int bx48=wMK->wEve.bx48;
117  int bx7=wMK->wEve.bx7;
118  if(spinDb->offsetBX48minusBX7(bx48,bx7)) {
119  printf("BAD bx7=%d bx48=%d del=%d\n",bx7,bx48,spinDb->offsetBX48minusBX7(bx48,bx7));
120  hA[0]->Fill("badBx48",1.);
121  return; // both counters must be in sync
122  }
123 
124  //remove events tagged as Zs
125  if(wMK->wEve.zTag) return;
126  hA[0]->Fill("noZ",1.);
127 
128 
129  hA[1]->Fill(bx48);
130  hA[2]->Fill(bx7);
131 
132  int bxStar48= spinDb->BXstarUsingBX48(bx48);
133  int bxStar7=spinDb->BXstarUsingBX7(bx7);
134  hA[3]->Fill(bxStar48);
135  hA[4]->Fill(bxStar7);
136 
137  int spin4=spinDb->spin4usingBX48(bx48);
138  hA[5]->Fill(bxStar7,spin4);
139 
140  float par_maxDsmThr=58;
141  float par_myET=25; // monitoring cut
142  if( wMK->wEve.l2bitRnd) { // lumi monitor BHT3-random
143  StMuEvent* muEve = wMK->mMuDstMaker->muDst()->event();
144  int max=0;
145  for (int m=0;m<300;m++) { // access L0-HT data
146  int val=muEve->emcTriggerDetector().highTower(m);
147  if(max<val) max=val;
148  }
149  // avoid too much energy - can be W-events (1/milion :)
150  if(max<par_maxDsmThr) { hA[6]->Fill(spin4); hA[0]->Fill("BG1",1.);}
151  return;
152  }
153 
154  if( wMK->wEve.l2bitET==0) return;
155  //..... it is guaranteed ..... L2W-ET>13 did fired ......
156 
157 
158  // search for Ws ............
159  for(unsigned int iv=0;iv<wMK->wEve.vertex.size();iv++) {
160  WeveVertex &V=wMK->wEve.vertex[iv];
161  for(unsigned int it=0;it<V.eleTrack.size();it++) {
162  WeveEleTrack &T=V.eleTrack[it];
163  if(T.pointTower.id==0) continue;
164 
165  /* Collect QCD background for lumi monitors */
166  float frac24=T.cluster.ET/(T.cl4x4.ET);
167  if(iv==0 && it==0 && frac24<wMK->par_clustFrac24) {
168  hA[31]->Fill(T.cluster.ET);
169  if( T.cluster.ET <20. ) { hA[7]->Fill(spin4); hA[0]->Fill("BG2",1.);}
170  }
171 
172  if(T.isMatch2Cl==false) continue;
173  assert(T.cluster.nTower>0); // internal logical error
174  assert(T.nearTotET>0); // internal logical error
175 
176  int iQ=0; // plus
177  float p_Q=T.prMuTrack->charge();
178  if( p_Q<0 ) iQ=1;// minus
179  float ET=T.cluster.ET;
180 
181 
182  //put final W cut here
183  bool isW= T.cluster.ET /T.nearTotET> wMK->par_nearTotEtFrac; // near cone
184  if(par_useNoEEMC)
185  isW=isW && T.sPtBalance_noEEMC>wMK->par_ptBalance; // awayET
186  else
187  isW=isW && T.sPtBalance>wMK->par_ptBalance; // awayET
188 
189  if(!isW) { // AL(QCD)
190  if(ET>15 &&ET<20 ) hA[16+iQ]->Fill(spin4);
191  continue;
192  }
193 
194  hA[0]->Fill("Wcut",1.);
195 
196  hA[30]->Fill(T.prMuTrack->eta());
197  if(fabs(T.primP.Eta()) > wMK->par_leptonEta) continue;
198  // allows further cuts on eta
199  if(T.prMuTrack->eta()<par_leptonEta1) continue;
200  if(T.prMuTrack->eta()>par_leptonEta2) continue;
201  hA[0]->Fill("eta",1.);
202 
203  //::::::::::::::::::::::::::::::::::::::::::::::::
204  //:::::accepted W events for x-section :::::::::::
205  //::::::::::::::::::::::::::::::::::::::::::::::::
206 
207  if(ET>par_myET) hA[0]->Fill("W25",1.);
208  float q2pt=T.prMuTrack->charge()/T.prMuTrack->pt();
209  if(ET>par_myET) hA[8]->Fill(q2pt);
210  hA[9]->Fill(ET,q2pt);
211 
212  // apply cut on reco charge
213  if( fabs(q2pt)< par_QPTlow) continue;
214  if(ET>par_myET) hA[0]->Fill("Qlow",1.);
215 
216  if(par_QPTlow>0) { // abaility to skip all Q/PT cuts
217  if( fabs(q2pt)< par_QPTlow) continue;
218  float highCut=par_QPThighA - (ET-par_QPThighET0)*par_QPThighB;
219  // printf("fff ET=%f q2pr=%f highCut=%f passed=%d\n",ET, q2pt,highCut,fabs(q2pt)<highCut);
220  if( ET>par_myET && ET<par_QPThighET1 && fabs(q2pt)>highCut) continue;
221  }
222 
223  if(ET>par_myET) {
224  hA[0]->Fill("Qhigh",1.);
225  if(p_Q>0) hA[0]->Fill("Q +",1.);
226  else hA[0]->Fill("Q -",1.);
227  }
228 
229 
230  hA[10+iQ]->Fill(ET);
231  if(ET>25 &&ET<50 ) hA[12+iQ]->Fill(spin4);
232  if(ET>32 &&ET<44 ) hA[14+iQ]->Fill(spin4);
233 
234  hA[18+iQ]->Fill(spin4,ET);
235 
236  } // end of loop over tracks
237  }// end of loop ove vertices
238 
239 }
240 
241 
242 // $Log: St2009pubSpinMaker.cxx,v $
243 // Revision 1.12 2010/05/01 01:31:45 balewski
244 // added W->JJ code & JES calibration
245 //
246 // Revision 1.11 2010/04/14 20:00:08 balewski
247 // added AL w/o endcap
248 //
249 // Revision 1.10 2010/03/22 17:18:32 balewski
250 // now the < > sign is right
251 //
252 // Revision 1.9 2010/03/22 16:11:42 balewski
253 // better computation of AL(QCD)
254 //
255 // Revision 1.8 2010/03/20 19:19:05 balewski
256 // added ability to drop Q/PT cut for spin analysis
257 //
258 // Revision 1.7 2010/03/15 17:05:46 balewski
259 // cleanup, used for W AL sort March 15, 2010
260 //
261 // Revision 1.6 2010/03/14 22:50:31 balewski
262 // *** empty log message ***
263 //
264 // Revision 1.5 2010/02/04 03:48:12 balewski
265 // add ET for lumi monitor
266 //
267 // Revision 1.4 2010/01/28 20:10:05 balewski
268 // added eta dependent spin sorting
269 //
270 // Revision 1.3 2010/01/28 03:42:55 balewski
271 // cleanup
272 //
273 // Revision 1.2 2010/01/27 22:12:25 balewski
274 // spin code matched to x-section code
275 //
276 // Revision 1.1 2009/11/23 23:00:18 balewski
277 // code moved spin-pool
278 //
279 // Revision 1.1 2009/11/23 21:11:18 balewski
280 // start
281 //
bool isBXfilledUsingInternalBX(int bx)
defined only for 2005 run by CAD , based on first 4 filled bunches in both rings. Note...
bool isValid()
dump spinDb content for current time stamp
Definition: StSpinDbMaker.h:53
StMuDst * muDst()
Definition: StMuDstMaker.h:425
Double_t pt() const
Returns pT at point of dca to primary vertex.
Definition: StMuTrack.h:256
int spin4usingBX48(int bx48)
8bit spin information
Short_t charge() const
Returns charge.
Definition: StMuTrack.h:255
int BX7offset()
bXing at STAR IP, [0,119]
Double_t eta() const
Returns pseudo rapidity at point of dca to primary vertex.
Definition: StMuTrack.h:257
void print(int level=0)
vs. STAR==yellow bXing
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
Definition: Stypes.h:40
int BX48offset()
bXing at STAR IP, [0,119]
int BXstarUsingBX7(int bx7)
bXing at STAR IP, [0,119]
spin sorting of W&#39;s
bool isPolDirLong()
Returns true if beams are transversely polarized, false otherwise.
Definition: StSpinDbMaker.h:55