StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSpaceChargeEbyEMaker.cxx
1 // //
3 // StSpaceChargeEbyEMaker performs event-by-event determination //
4 // of the space charge correction for tracks, and sets it for //
5 // the next event. //
6 // //
8 
9 #include "TROOT.h"
10 #include "TSystem.h"
11 #include "StSpaceChargeEbyEMaker.h"
12 #include "StDbUtilities/StMagUtilities.h"
13 #include "StEvent/StEventTypes.h"
14 #include "StMessMgr.h"
15 #include "StTpcDb/StTpcDb.h"
16 #include "StTpcDb/StTpcDbMaker.h"
17 #include "St_db_Maker/St_db_Maker.h"
18 #include "StEvent/StDcaGeometry.h"
19 #include "StEvent/StBTofCollection.h"
20 #include "StEvent/StBTofHit.h"
21 #include "StEvent/StBTofHeader.h"
22 #include "StEvent/StBTofPidTraits.h"
23 #include "StEvent/StTrackPidTraits.h"
24 #include "StDetectorDbMaker/St_trigDetSumsC.h"
25 #include "StDetectorDbMaker/St_tpcGridLeakC.h"
26 #include "StDetectorDbMaker/St_spaceChargeCorC.h"
27 #include "StEmcUtil/geometry/StEmcGeom.h"
28 #include "StEmcUtil/projection/StEmcPosition.h"
29 
30 #include "TUnixTime.h"
31 #include "TFile.h"
32 #include "TH2.h"
33 #include "TH3.h"
34 #include "TF1.h"
35 #include "TNtuple.h"
36 #include "TPad.h"
37 
38 // Histogram ranges:
39 const int SCN1 = 400;
40 const int SCN2 = 100;
41 const float SCL = -0.015;
42 const float SCH = 0.185;
43 const float SCB = (SCH-SCL)/SCN1;
44 const float SCX = SCB/TMath::Sqrt(TMath::TwoPi());
45 
46 const int DCN = 125;
47 const float DCL = -2.5;
48 const float DCH = 2.5;
49 
50 const int PHN = 72;
51 const float PI2 = TMath::TwoPi();
52 
53 const int EVN = 1024;
54 
55 const int ZN = 60;
56 const float ZL = -150.;
57 const float ZH = 150.;
58 
59 const int GN = 150;
60 const float GL = -0.3;
61 const float GH = 1.2;
62 const int GZN = 17;
63 const float GZL = -200.;
64 const float GZH = 225.;
65 
66 const float ZGGRID = 208.707;
67 
68 static TF1 ga1("ga1","[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))");
69 static TF1 ln1("ln1","[0]+((208.707-abs(x))*[1]/100.0)",-150.,150.);
70 
71 const unsigned int MAXVTXCANDIDATES = 256;
72 
73 ClassImp(StSpaceChargeEbyEMaker)
74 
75 //_____________________________________________________________________________
77  event(0),runinfo(0),
78  Calibmode(kFALSE), PrePassmode(kFALSE), PrePassdone(kFALSE),
79  QAmode(kFALSE), TrackInfomode(0), Asymmode(kFALSE),
80  doNtuple(kFALSE), doReset(kTRUE), doGaps(kFALSE), doSecGaps(kFALSE),
81  inGapRow(0),
82  vtxVpdAgree(5.0), vtxPCTs(0), vtxEmcMatch(1), vtxTofMatch(0),
83  vtxMinTrks(5), minTpcHits(25),
84  reqEmcMatch(kFALSE), reqTofMatch(kFALSE), reqEmcOrTofMatch(kTRUE),
85  m_ExB(0), SCcorrection(0), GLcorrection(0), SCEWRatio(0),
86  scehist(0), timehist(0), myhist(0), myhistN(0), myhistP(0),
87  myhistE(0), myhistW(0), dczhist(0), dcehist(0), dcphist(0),
88  dcahist(0), dcahistN(0), dcahistP(0), dcahistE(0), dcahistW(0),
89  gapZhist(0), gapZhistneg(0), gapZhistpos(0), cutshist(0), ntup(0) {
90 
91  MINTRACKS=6000;
92  //SCALER_ERROR = 0.0006; // by eye from hist: SCvsZDCEpW.gif (liberal)
93  SCALER_ERROR = 0.0007; // by RMS from hist: SCvsZDCX.gif (liberal)
94 
95  // MAXDIFFE is maximum different in sc from last ebye sc
96  MAXDIFFE = SCALER_ERROR;
97  // MAXDIFFA is maximum different in sc from last scaler sc
98  MAXDIFFA = 2*SCALER_ERROR; // should be about equal to SCALER_ERROR, no?
99  // Present uncetainties with scalers demands greater tolerance
100 
101  // initializations needed at the start
102  runid = 0;
103  memset(evts,0,SCHN*sizeof(int));
104  memset(times,0,SCHN*sizeof(int));
105  memset(evtstbin,0,SCHN*sizeof(float));
106  evtsnow = 0;
107  firstEvent = -1;
108 
109  SetMode(0); // default is mode 0 (no QA, no PrePass)
110  //DoQAmode(); // For testing
111 
112  schist = new TH1F("SpCh","Space Charge" ,SCN1,SCL,SCH);
113  schistE = new TH1F("SpChE","Space Charge East",SCN1,SCL,SCH);
114  schistW = new TH1F("SpChW","Space Charge West",SCN1,SCL,SCH);
115  schist ->SetDirectory(0);
116  schistE->SetDirectory(0);
117  schistW->SetDirectory(0);
118  for (int i=0;i<SCHN;i++){
119  schists[i] = new TH1F(Form("SpCh%d" ,i),"Space Charge" ,SCN1,SCL,SCH);
120  schistsE[i] = new TH1F(Form("SpChE%d",i),"Space Charge East",SCN1,SCL,SCH);
121  schistsW[i] = new TH1F(Form("SpChW%d",i),"Space Charge West",SCN1,SCL,SCH);
122  schists[i] ->SetDirectory(0);
123  schistsE[i]->SetDirectory(0);
124  schistsW[i]->SetDirectory(0);
125  }
126  for (int i=0; i<24; i++) {
127  schistS[i] = 0;
128  gapZhistS[i] = 0;
129  gapZhistnegS[i] = 0;
130  gapZhistposS[i] = 0;
131  }
132 
133  // other initializations for safety (will be set later anyhow)
134  evt = 0;
135  curhist = 0;
136  lasttime = 0;
137  memset(scS,0,24*sizeof(float));
138  memset(escS,0,24*sizeof(float));
139  sc = 0; esc = 0;
140  scE = 0; escE = 0;
141  scW = 0; escW = 0;
142  lastsc = 0;
143  lastEWRatio = 0;
144  oldevt = 0;
145  did_auto = kTRUE;
146  memset(ntrks ,0,SCHN*sizeof(float));
147  memset(ntrksE,0,SCHN*sizeof(float));
148  memset(ntrksW,0,SCHN*sizeof(float));
149  memset(ntrksS,0,24*sizeof(float));
150  gapZfitslope = 0; gapZfitintercept = 0; gapZdivslope = 0;
151  egapZfitslope = 0; egapZfitintercept = 0; egapZdivslope = 0;
152  gapZfitslopeneg = 0; gapZfitinterceptneg = 0; gapZdivslopeneg = 0;
153  gapZfitslopepos = 0; gapZfitinterceptpos = 0; gapZdivslopepos = 0;
154  gapZfitslopeeast = 0; gapZfitintercepteast = 0; gapZdivslopeeast = 0;
155  gapZfitslopewest = 0; gapZfitinterceptwest = 0; gapZdivslopewest = 0;
156  memset(gapZfitslopeS,0,24*sizeof(float));
157  memset(gapZfitinterceptS,0,24*sizeof(float));
158  memset(gapZdivslopeS,0,24*sizeof(float));
159  memset(gapZfitslopenegS,0,24*sizeof(float));
160  memset(gapZfitinterceptnegS,0,24*sizeof(float));
161  memset(gapZdivslopenegS,0,24*sizeof(float));
162  memset(gapZfitslopeposS,0,24*sizeof(float));
163  memset(gapZfitinterceptposS,0,24*sizeof(float));
164  memset(gapZdivslopeposS,0,24*sizeof(float));
165 }
166 //_____________________________________________________________________________
167 StSpaceChargeEbyEMaker::~StSpaceChargeEbyEMaker() {
168  delete schist;
169  delete schistE;
170  delete schistW;
171  for (int i=0;i<SCHN;i++) {
172  delete schists[i];
173  delete schistsE[i];
174  delete schistsW[i];
175  }
176 }
177 //_____________________________________________________________________________
179 
180  if (evt > 0) {
181  if (PrePassmode) {
182  if (PrePassdone) WriteTableToFile();
183  else gMessMgr->Warning("StSpaceChargeEbyEMaker: NO SC => NO OUTPUT");
184  }
185  if ((!Calibmode)&&(!PrePassdone)) EvalCalib();
186  WriteQAHists();
187  } else {
188  gMessMgr->Warning("StSpaceChargeEbyEMaker: NO EVENTS => NO OUTPUT");
189  }
190 
191  return kStOk;
192 }
193 //_____________________________________________________________________________
194 Int_t StSpaceChargeEbyEMaker::Init() {
195 
196  // Use mode to set switches:
197  // Set mode in BFC chain by attribute goptSCE100XXX, where
198  // XXX is the mode number, e.g. goptSCE100050 sets mode 50
199  Int_t attrMode = IAttr(".gopt.sce");
200  attrMode = (attrMode ? attrMode%1000 : GetMode());
201  gMessMgr->Info() << "StSpaceChargeEbyEMaker mode: " << attrMode << endm;
202  switch (attrMode) {
203  case (1) : DoQAmode(); break;
204  case (2) : DoPrePassmode(); break;
205  case (3) : DoPrePassmode(); DoQAmode(); break;
206  case (4) : DoCalib(); break;
207  case (5) : DoCalib(); DoAsym(); break;
208  case (6) : DoCalib(); DoSecGaps(); break;
209  case (7) : DoCalib(); DoAsym(); DoSecGaps(); break;
210  case (10): DoNtuple(); break;
211  case (11): DoNtuple(); DontReset(); break;
212  case (12): DoNtuple(); DoQAmode(); break;
213  case (13): DoNtuple(); DontReset(); DoQAmode(); break;
214  case (20): DoNtuple(); DontReset(); DoQAmode(); DoSecGaps(); break;
215  case (50): DoTrackInfo(); break; // filter out pile-up tracks
216  case (51): DoTrackInfo(2); break; // include pile-up tracks
217  case (52): DoTrackInfo(3); break; // include pile-up tracks of any length
218  default : {}
219  }
220 
221  if (Calibmode) doReset = kFALSE;
222 
223  evt=0;
224  oldevt=1;
225  lastsc=0.;
226  lastEWRatio=0.;
227  curhist=0;
228  lasttime=0;
229  did_auto=kTRUE;
230  InitQAHists();
231  if (QAmode) gMessMgr->Info("StSpaceChargeEbyEMaker: Initializing");
232  if (TrackInfomode) gMessMgr->Info("StSpaceChargeEbyEMaker: Track Info mode");
233  return StMaker::Init();
234 }
235 //_____________________________________________________________________________
237 
238  if (QAmode) cutshist->Fill(0);
239 
240  // On very first event, determine first event timestamp and
241  // set default parameters, unless in Calibmode
242  if ((!Calibmode) && (tabname.Length() == 0)) SetTableName();
243 
244  // Get instance of StMagUtilities
245  m_ExB = StMagUtilities::Instance();
246  if (!m_ExB) {
247 #ifdef __NEW_MagUtilities__
248  m_ExB = new StMagUtilities(gStTpcDb,(kSpaceChargeR2 | kGridLeak));
249 #else /* ! __NEW_MagUtilities__ */
250  TDataSet *RunLog = GetDataBase("RunLog/MagFactor");
251  if (!RunLog) gMessMgr->Warning("StSpaceChargeEbyEMaker: No RunLog/MagFactor found.");
252  m_ExB = new StMagUtilities(gStTpcDb,RunLog,(kSpaceChargeR2 | kGridLeak));
253 #endif /* __NEW_MagUtilities__ */
254  }
255  m_ExB->UndoDistortion(0,0,0); // initialize for this event in case it wasn't used
256  lastsc = m_ExB->CurrentSpaceChargeR2();
257  lastEWRatio = m_ExB->CurrentSpaceChargeEWRatio();
258 
259  // Get StEvent and related info, determine if things are OK
260  event = static_cast<StEvent*>(GetInputDS("StEvent"));
261  if (!event) {
262  gMessMgr->Warning("StSpaceChargeEbyEMaker: no StEvent; skipping event.");
263  return kStWarn;
264  }
265  if (QAmode) cutshist->Fill(1);
266  if (firstEvent<0) firstEvent = event->id();
267  // Get runinfo, determine if the magnetic field is nonzero
268  // EbyE maker not currently able to handle zero B field
269  runinfo = event->runInfo();
270  if ((!runinfo) || (runinfo->magneticField() == 0)) {
271  gMessMgr->Error("StSpaceChargeEbyEMaker: cannot run due to zero or unknown mag field!");
272  // Look for any errant zero field SC values as a warning
273  // for processing even without EbyE
274  if ((lastsc != 0) && (runinfo) && (runinfo->magneticField() == 0))
275  gMessMgr->Warning() << "BE AWARE THAT A NONZERO VALUE OF SPACECHARGE\n"
276  << " WAS RETURNED BY DB CALL!\n (could be a local DB file or"
277  << " in the actual database)" << endm;
278  return kStFatal;
279  }
280  if (QAmode) cutshist->Fill(2);
281 
282  // Select highest ranked vertex(vertices) + some quality cuts
283  StPrimaryVertex* pvtx = 0;
284  unsigned int vtxCandidates[MAXVTXCANDIDATES];
285  unsigned int numVtxCandidates = 0;
286  unsigned int totVertices = event->numberOfPrimaryVertices();
287  if (TrackInfomode>1) numVtxCandidates=1;
288  else {
289  const StBTofCollection* btofColl = event->btofCollection();
290  const StBTofHeader* btofHeader = (btofColl ? btofColl->tofHeader() : 0);
291  float vpd_zvertex = (btofHeader ? btofHeader->vpdVz() : -999);
292  for (unsigned int vtxIdx = 0; vtxIdx < totVertices; vtxIdx++) {
293  pvtx = event->primaryVertex(vtxIdx);
294  if (QAmode) cutshist->Fill(3);
295  if (! (IAttr("EastOff") || IAttr("WestOff"))) {
296  // vertex ranking & ordering break for East/West off
297  StVertexFinderId vtxFindID = pvtx->vertexFinderId();
298  float min_rank = -1e6;
299  switch (vtxFindID) {
300  case minuitVertexFinder : min_rank = -5; break;
301  case ppvVertexFinder :
302  case ppvNoCtbVertexFinder : min_rank = 0; break;
303  default : break;
304  }
305  // only one chance for MinuitVF
306  if (vtxFindID == minuitVertexFinder) totVertices = 1;
307  // vertices are rank ordered, so once it fails, we're done
308  if (pvtx->ranking() < min_rank) break;
309  }
310  if (QAmode) cutshist->Fill(4);
311  if (pvtx->numberOfDaughters() < vtxMinTrks) continue;
312  if (QAmode) cutshist->Fill(5);
313  if (pvtx->numMatchesWithBEMC() < vtxEmcMatch) continue;
314  if (QAmode) cutshist->Fill(6);
315  if (pvtx->numMatchesWithBTOF() < vtxTofMatch) continue;
316  if (QAmode) cutshist->Fill(7);
317  if (pvtx->numPostXTracks() > vtxPCTs) continue;
318  if (QAmode) cutshist->Fill(8);
319  if (vtxVpdAgree > 0 && // set vtxVpdAgree negative to skip this cut
320  TMath::Abs(pvtx->position().z() - vpd_zvertex) > vtxVpdAgree) continue;
321  if (QAmode) cutshist->Fill(9);
322  vtxCandidates[numVtxCandidates] = vtxIdx;
323  numVtxCandidates++;
324  if (numVtxCandidates == MAXVTXCANDIDATES) break;
325  }
326  }
327  if (!numVtxCandidates) return kStOk;
328  if (QAmode) cutshist->Fill(10);
329 
330  StSPtrVecTrackNode& theNodes = event->trackNodes();
331  unsigned int nnodes = theNodes.size();
332  if (!nnodes) return kStOk;
333  if (QAmode) cutshist->Fill(11);
334 
335  // Store and setup event-wise info
336  evt++;
337  int thistime = event->time();
338  if (lasttime) {
339  timehist->SetBinContent(evt,thistime-lasttime);
340  } else {
341  runid = event->runId();
342  }
343  if (doReset) {
344  if (evt>1) curhist = imodHN(curhist+1);
345  schists[curhist] ->Reset();
346  schistsE[curhist]->Reset();
347  schistsW[curhist]->Reset();
348  if (doGaps) {
349  gapZhist->Reset();
350  gapZhistpos->Reset();
351  gapZhistneg->Reset();
352  if (doSecGaps) for (int i=0; i<24; i++) {
353  schistS[i]->Reset();
354  gapZhistS[i]->Reset();
355  gapZhistposS[i]->Reset();
356  gapZhistnegS[i]->Reset();
357  }
358  }
359  } else {
360  // Do not reset ntuple in calib mode
361  if (doNtuple && !Calibmode) ntup->Reset();
362  }
363 
364  // Keep time and event number
365  times[curhist] = thistime;
366  evts[curhist]=evt;
367 
368  // Keep calibrations used
369  if (!SCcorrection) {
370  SCcorrection = new TNamed("SCcorrection",
371  (St_spaceChargeCorR2C::instance()->getSpaceChargeString()).Data());
372  SCEWRatio = new TNamed("SCEWRatio",Form("%f",
373  St_spaceChargeCorR2C::instance()->getEWRatio()));
374  GLcorrection = new TNamed("GLcorrection",Form("%f",
375  St_tpcGridLeakC::instance()->MiddlGLStrength()));
376  gMessMgr->Info() << "Using the following corrections:" << endm;
377  gMessMgr->Info() << "sc = " << SCcorrection->GetTitle() << endm;
378  gMessMgr->Info() << "E/W= " << SCEWRatio->GetTitle() << endm;
379  gMessMgr->Info() << "GL = " << GLcorrection->GetTitle() << endm;
380  }
381 
382  // Keep track of # of events in the same time bin
383  if (thistime == lasttime) evtsnow++;
384  else evtsnow = 1;
385  evtstbin[curhist] = evtsnow;
386 
387  if (QAmode) {
388  gMessMgr->Info()
389  << "used (for this event) SpaceCharge = "
390  << lastsc << " (" << thistime << ")" << endm;
391  gMessMgr->Info()
392  << "zdc west+east = "
393  << runinfo->zdcWestRate()+runinfo->zdcEastRate() << endm;
394  gMessMgr->Info()
395  << "zdc coincidence = "
396  << runinfo->zdcCoincidenceRate() << endm;
397  }
398 
399  // Fill the StEvent information for the SpaceCharge used in this event
400  runinfo->setSpaceCharge(lastsc);
401  runinfo->setSpaceChargeCorrectionMode(m_ExB->GetSpaceChargeMode());
402  if (!inGapRow) {
403  if (runinfo->runId() > 10000000) inGapRow = 13; // Run 9+
404  else if (runinfo->runId() > 0) inGapRow = 12; // Run exists
405  // else undefined
406  }
407 
408  // Track loop
409  unsigned int i,j,k,v;
410 
411  // Prepare for EMC match
412  StEmcDetector* bemcDet = 0;
413  Double_t emcRadius = 0;
414  static StEmcPosition* emcPosition = 0;
415  static StEmcGeom* emcGeom = 0;
416  if (reqEmcMatch || reqEmcOrTofMatch) {
417  bemcDet = event->emcCollection()->detector(kBarrelEmcTowerId);
418  if (!emcPosition) emcPosition = new StEmcPosition();
419  if (!emcGeom) emcGeom = StEmcGeom::instance("bemc");
420  emcRadius = emcGeom->Radius() + 30; // use exit radius, 30cm beyond face
421  }
422 
423  St_tpcPadConfigC* pads = St_tpcPadConfigC::instance();
424  StThreeVectorD vtxPos(0,0,0),vtxPosErr(0,0,0);
425  for (v=0; v<numVtxCandidates; v++) {
426  if (TrackInfomode<2) {
427  pvtx = event->primaryVertex(vtxCandidates[v]);
428  if (! pvtx) continue;
429  vtxPos = pvtx->position();
430  vtxPosErr = pvtx->positionError();
431  }
432 
433  for (i=0; i<nnodes; i++) {
434  for (j=0; j<theNodes[i]->entries(global); j++) {
435  if (QAmode) cutshist->Fill(16);
436  StTrack* tri = theNodes[i]->track(global,j);
437  if (!tri) continue;
438  if (QAmode) cutshist->Fill(17);
439 
440  const StTrackTopologyMap& map = tri->topologyMap();
441  //if (! map.trackTpcOnly()) continue;
442  if (! map.hasHitInDetector(kTpcId)) continue;
443  if (QAmode) cutshist->Fill(18);
444  // Multiple silicon hits destroy sDCA <-> SpaceCharge correlation,
445  // and single hit in SVT is unreliable. Only good config is NO SVT!
446  // Updated for HFT era to exclude anything with PXL or IST hits
447  if ((map.hasHitInDetector(kSvtId) ||
448  map.hasHitInDetector(kPxlId) ||
449  map.hasHitInDetector(kIstId)) && TrackInfomode < 1) continue;
450  if (QAmode) cutshist->Fill(19);
451  if (map.numberOfHits(kTpcId) < minTpcHits) continue;
452  if (QAmode) cutshist->Fill(20);
453 
454  // *** TOF MATCHING ***
455  Bool_t tofMatch = kFALSE;
456  if (reqTofMatch || reqEmcOrTofMatch) {
457  const StPtrVecTrackPidTraits& theTofPidTraits = tri->pidTraits(kTofId);
458  if (theTofPidTraits.size()) {
459  if (QAmode) cutshist->Fill(21);
460  StTrackPidTraits* theSelectedTrait = theTofPidTraits[theTofPidTraits.size()-1];
461  if (theSelectedTrait) {
462  if (QAmode) cutshist->Fill(22);
463  StBTofPidTraits *pidTof = dynamic_cast<StBTofPidTraits *>(theSelectedTrait);
464  if (pidTof) {
465  if (QAmode) cutshist->Fill(23);
466  int Mflag=pidTof->matchFlag();
467  // 0: no matching
468  // 1: 1-1 matching
469  // 2: 1-2 matching, pick up the one with higher ToT vaule (<25ns)
470  tofMatch = (Mflag > 0);
471  }
472  }
473  }
474  // Don't require match if reqEmcOrTofMatch
475  if (reqTofMatch && !tofMatch) continue;
476  }
477  if (QAmode) cutshist->Fill(24);
478 
479  // *** EMC MATCHING ***
480  Bool_t emcMatch = kFALSE;
481  if (reqEmcMatch || (reqEmcOrTofMatch&&(!tofMatch))) {
482  Double_t mEmcThresh = 0.15;
483  Double_t energyBEMC = -100.0;
484  UInt_t tower_eta,tower_mod = 0;
485  Int_t tower_sub = 0;
486  StThreeVectorD emcTrkMomentum,emcTrkPosition;
487  if ((emcPosition->trackOnEmc(&emcTrkPosition,&emcTrkMomentum,
488  tri,runinfo->magneticField()/10.,emcRadius))) {
489  if (QAmode) cutshist->Fill(25);
490  Float_t emcEta = emcTrkPosition.pseudoRapidity();
491  Float_t emcPhi = emcTrkPosition.phi();
492  Int_t m = 0 ,e = 0,s = 0,id = 0;
493  emcGeom->getBin(emcPhi,emcEta,m,e,s);
494  if (emcGeom->getId(m,e,s,id) == 0) {
495  tower_mod = m;
496  tower_eta = e;
497  tower_sub = s;
498  }
499  if (tower_mod >= 1 && tower_mod <= 120) {
500  if (QAmode) cutshist->Fill(26);
501  if (event->emcCollection()) {
502  StEmcModule* emcMod = bemcDet->module(tower_mod);
503  StSPtrVecEmcRawHit& emcHits = emcMod->hits();
504  for (UInt_t emcHit=0; emcHit<emcHits.size(); emcHit++) {
505  if ((emcHits[emcHit]) && (emcHits[emcHit]->eta() == tower_eta)
506  && (emcHits[emcHit]->sub() == tower_sub)) {
507  energyBEMC = emcHits[emcHit]->energy();
508  break; // only one hit
509  }
510  }
511  emcMatch = (energyBEMC >= mEmcThresh);
512  }
513  }
514  }
515  // Require match for both reqEmcMatch or reqEmcOrTofMatch if here
516  if (!emcMatch) continue;
517  }
518  if (QAmode) cutshist->Fill(27);
519 
520  StTrackGeometry* triGeom = tri->geometry();
521 
522  StThreeVectorF xvec = triGeom->origin();
523  if (!(xvec.x() || xvec.y() || xvec.z())) continue;
524  if (QAmode) cutshist->Fill(28);
525  StThreeVectorF pvec = triGeom->momentum();
526  if (!(pvec.x() || pvec.y())) continue;
527  if (QAmode) cutshist->Fill(29);
528 
529  float oldPt = pvec.perp();
530  if (oldPt < 0.0001) continue;
531  if (QAmode) cutshist->Fill(30);
532 
533  int e_or_w = 0; // east is -1, west is +1
534  if (pvec.z() * xvec.z() > 0) e_or_w = ( (xvec.z() > 0) ? 1 : -1 );
535 
536  StPhysicalHelixD hh = triGeom->helix();
537 
538  Float_t eta=pvec.pseudoRapidity();
539  Float_t phi=0;
540  Double_t pathlen = 0;
541  //Float_t DCA=hh.geometricSignedDistance(0,0); // for testing only
542  Float_t DCA3=-999;
543  Float_t DCA2=-999;
544  Double_t DCAerr = 0.;
545  StDcaGeometry* triDcaGeom = ((StGlobalTrack*) tri)->dcaGeometry();
546  if (triDcaGeom) {
547  StPhysicalHelixD dcahh = triDcaGeom->helix();
548  DCA3 = dcahh.distance(vtxPos,kFALSE);
549  DCA2 = dcahh.geometricSignedDistance(vtxPos.x(),vtxPos.y());
550  // helix() gets the sign of DCA2, thelix() gets the error
551  THelixTrack thelix = triDcaGeom->thelix();
552  thelix.Dca(vtxPos.x(),vtxPos.y(),&DCAerr);
553  phi = TMath::ATan2(dcahh.cy(pathlen),dcahh.cx(pathlen));
554  if (TrackInfomode>1) {
555  vtxPos.setZ(dcahh.z(pathlen));
556  }
557  } else {
558  DCA3 = hh.distance(vtxPos,kFALSE);
559  DCA2 = hh.geometricSignedDistance(vtxPos.x(),vtxPos.y());
560  pathlen = hh.pathLength(vtxPos.x(),vtxPos.y());
561  phi = TMath::ATan2(hh.cy(pathlen),hh.cx(pathlen));
562  if (TrackInfomode>1) {
563  vtxPos.setZ(hh.z(pathlen));
564  vtxPosErr.setY(1); // flag the non-DcaGeom tracks
565  }
566  }
567  if (TrackInfomode>1) {
568  if (TMath::Abs(DCA2) > 2) continue; // cut out tracks not near (0,0)
569  } else {
570  if (DCA3 > 4) continue; // cut out pileup tracks!
571  }
572  if (QAmode) cutshist->Fill(31);
573  Int_t ch = (int) triGeom->charge();
574 
575  Int_t PCT = 0;
576  Int_t sec,sector = -1;
577  Int_t prow = 0;
578  Int_t prowmask = 0x0;
579  Double_t rs[128];
580  Double_t rerrors[128];
581  Double_t rphierrors[128];
582  memset(rerrors,0,64*sizeof(Float_t));
583  memset(rphierrors,0,64*sizeof(Float_t));
584  StPtrVecHit& hits = tri->detectorInfo()->hits();
585  for (k=0;k<hits.size();k++) {
586  StHit* hit = hits[k];
587  if (hit->detector() == kTpcId) {
588  sec = ((StTpcHit*) hit)->sector();
589  prow = ((StTpcHit*) hit)->padrow();
590  if ((hit->position().z() > 1 && sec > 12) ||
591  (hit->position().z() <-1 && sec < 13)) PCT++;
592  int lastInner = pads->innerPadRows(sec);
593  if (prow >= lastInner-2 && prow <= lastInner+3) {
594  sector = ( sec == sector || sector < 0 ? sec : 0 ); // 0 if crossing sectors
595  // Require 4 hits: one on prow closest to gap for both inner and outer,
596  // and another on either of next two prows away from gap for both inner and outer
597  // (exception for inGapRow==12 to require prow 11 & 12 for inner)
598  int shifter = prow-(lastInner-1);
599  if (shifter<0) shifter += (inGapRow==12 ? 2 : 1);
600  else if (shifter>3) shifter--;
601  prowmask |= (0x1 << shifter);
602  }
603  }
604  rs[k] = hit->position().perp();
605  Float_t herr = hit->positionError().perp();
606  rerrors[k] = herr;
607  rphierrors[k] = herr;
608  }
609  if (PCT && TrackInfomode<2) continue; // Track has post-crossing hits
610  bool good4gap = (sector > 0) && (prowmask == 0xF);
611  if (QAmode) cutshist->Fill(32);
612 
613  if (TrackInfomode) {
614  LOG_INFO << Form("GOODTRACK %d %d %6.2f %9.4f %8.3f %8.4f %8.4f %6.4f %6.4f %d",
615  runid,event->id(),vtxPos.z(),ch/oldPt,eta,phi,DCA2,DCAerr,
616  vtxPosErr.perp(),hits.size()) << endm;
617  continue;
618  }
619 
620  Float_t space = 10000.;
621  Int_t predictFailed = m_ExB->PredictSpaceChargeDistortion(hits.size(),ch,oldPt,vtxPos.z(),
622  eta,phi,DCA2,rs,rerrors,rphierrors,space);
623  if (predictFailed) {
624  if (QAmode) cutshist->Fill(40+predictFailed);
625  continue;
626  }
627  if (QAmode) cutshist->Fill(33);
628 
629  TH1F** aschists = (e_or_w > 0 ? schistsW :
630  (e_or_w < 0 ? schistsE : 0));
631  Bool_t doSecSCs = (doSecGaps && sector>0);
632  sector--;
633 
634  Double_t spaceErr = TMath::Abs(space*DCAerr/DCA2);
635  Float_t spaceEW = space + lastsc * (e_or_w < 0 ? lastEWRatio : 1.0);
636  space += lastsc; // Assumes additive linearity of space charge!
637  if (spaceErr > 0) {
638  // Fill SpaceCharge accounting for prediction error
639  ga1.SetParameters(SCX/spaceErr,space,spaceErr);
640  for (Int_t si=1;si<=SCN1;si++) {
641  Double_t sx = schists[curhist]->GetBinCenter(si);
642  if (TMath::Abs((sx-space)/spaceErr) < 4.5) {
643  // only within +/-4.5 sigma
644  schists[curhist]->Fill(sx,ga1.Eval(sx));
645  }
646  }
647  // Now for east/west
648  if (aschists || doSecSCs) {
649  ga1.SetParameters(SCX/spaceErr,spaceEW,spaceErr);
650  for (Int_t si=1;si<=SCN1;si++) {
651  Double_t sx = schists[curhist]->GetBinCenter(si);
652  if (TMath::Abs((sx-spaceEW)/spaceErr) < 4.5) {
653  // only within +/-4.5 sigma
654  if (aschists) aschists[curhist]->Fill(sx,ga1.Eval(sx));
655  if (doSecSCs) schistS [sector ]->Fill(sx,ga1.Eval(sx));
656  }
657  }
658  }
659  } else {
660  schists[curhist]->Fill(space);
661  if (aschists) aschists[curhist]->Fill(spaceEW);
662  if (doSecSCs) schistS [sector ]->Fill(spaceEW);
663  }
664  FillQAHists(DCA2,space,spaceEW,ch,hh,e_or_w);
665 
666 
667  if (doGaps && good4gap &&
668  (e_or_w!=0) && (TMath::Abs(ch)==1) && (oldPt>0.3))
669  FillGapHists(tri,hh,e_or_w,ch);
670 
671  } // loop over j tracks
672  } // loop over i Nodes
673  } // loop over v vertices
674  if (QAmode) cutshist->Fill(9);
675 
676 
677  ntrks[curhist] = schists[curhist] ->Integral();
678  ntrksE[curhist] = schistsE[curhist]->Integral();
679  ntrksW[curhist] = schistsW[curhist]->Integral();
680  if (doSecGaps) for (i=0; i<24; i++) ntrksS[i] = schistS[i]->Integral();
681 
682  // Wrap it up and make a decision
683  int result = DecideSpaceCharge(thistime);
684 
685  if (doGaps) DetermineGaps();
686  if (doNtuple) {
687  static float X[123];
688  static float ntent = 0.0;
689  static float nttrk = 0.0;
690 
691  if (ntent == 0.0) memset(X,0,51*sizeof(float));
692  ntent++; // # entries since last reset, including this one
693  float last_nttrk = nttrk;
694  nttrk = ntrks[curhist]; // # tracks since last reset, including these
695  float s0 = ( nttrk ? last_nttrk / nttrk : 0 );
696  float s1 = 1.0 - s0; // fraction of tracks from current event
697 
698  if (QAmode) {
699  gMessMgr->Info() << "reset = " << doReset << endm;
700  gMessMgr->Info() << "nevts = " << ntent << endm;
701  gMessMgr->Info() << "ntrks = " << nttrk << endm;
702  if (doSecGaps) for (i=0; i<24; i++) {
703  gMessMgr->Info() << "ntrks(" << i+1 << ") = " << ntrksS[i] << endm;
704  }
705  }
706 
707  float ee;
708  int fbin = evt + 1 - ((int) ntent);
709 
710  X[0] = sc;
711  X[1] = FindPeak(dcehist->ProjectionY("_dy",fbin,evt),ee);
712  X[2] = s0*X[2] + s1*runinfo->zdcCoincidenceRate();
713  X[3] = s0*X[3] + s1*runinfo->zdcWestRate();
714  X[4] = s0*X[4] + s1*runinfo->zdcEastRate();
715  X[5] = s0*X[5] + s1*runinfo->bbcCoincidenceRate();
716  X[6] = s0*X[6] + s1*runinfo->bbcWestRate();
717  X[7] = s0*X[7] + s1*runinfo->bbcEastRate();
718  X[8] = s0*X[8] + s1*runinfo->bbcBlueBackgroundRate();
719  X[9] = s0*X[9] + s1*runinfo->bbcYellowBackgroundRate();
720  X[10] = s0*X[10] + s1*runinfo->initialBeamIntensity(blue); // west-bound
721  X[11] = s0*X[11] + s1*runinfo->initialBeamIntensity(yellow); // east-bound
722  X[12] = runinfo->beamFillNumber(blue);
723  X[13] = runinfo->magneticField();
724  X[14] = event->runId();
725  X[15] = firstEvent;
726  if ((QAmode) && (evt <= EVN)) {
727  X[16] = FindPeak(dcahistN->ProjectionZ("_dnz",fbin,evt,1,PHN),ee);
728  X[17] = FindPeak(dcahistP->ProjectionZ("_dpz",fbin,evt,1,PHN),ee);
729  X[18] = FindPeak(dcahistE->ProjectionZ("_dez",fbin,evt,1,PHN),ee);
730  X[19] = FindPeak(dcahistW->ProjectionZ("_dwz",fbin,evt,1,PHN),ee);
731  }
732  X[20] = gapZfitslope;
733  X[21] = gapZfitintercept;
734  X[22] = gapZdivslope;
735  X[23] = gapZfitslopeneg;
736  X[24] = gapZfitinterceptneg;
737  X[25] = gapZdivslopeneg;
738  X[26] = gapZfitslopepos;
739  X[27] = gapZfitinterceptpos;
740  X[28] = gapZdivslopepos;
741  X[29] = gapZfitslopeeast;
742  X[30] = gapZfitintercepteast;
743  X[31] = gapZdivslopeeast;
744  X[32] = gapZfitslopewest;
745  X[33] = gapZfitinterceptwest;
746  X[34] = gapZdivslopewest;
747  X[35] = s0*X[35] + s1*runinfo->spaceCharge();
748  X[36] = s0*X[36] + s1*((float) (runinfo->spaceChargeCorrectionMode()));
749  X[37] = s0*X[37] + s1*St_tpcGridLeakC::instance()->MiddlGLStrength();
750  X[38] = s0*X[38] + s1*St_trigDetSumsC::Nc(runinfo->zdcCoincidenceRate(),
751  runinfo->zdcEastRate(),runinfo->zdcWestRate());
752  X[39] = s0*X[39] + s1*St_trigDetSumsC::Nc(runinfo->bbcCoincidenceRate(),
753  runinfo->bbcEastRate(),runinfo->bbcWestRate());
754 
755  // VPD data:
756  // StRunInfo's backgroundRate is filled in StEventMaker from 'mult' of trigDetSums
757  // trigDetSums fills 'mult' from rich scaler rs11 in the DAQ stream
758  // (but in the offline database, it is migrated from rs16)
759  // rs11 stores VPD coincidence rate as of 2007-12-19
760  // rs16 stores MTD rate as of Run 9
761  // VPD east and west are rs8 and 9, and are not stored in StEvent
762  X[40] = s0*X[40] + s1*runinfo->backgroundRate();
763  X[41] = s0*X[41] + s1*St_trigDetSumsC::instance()->getPVPDWest();
764  X[42] = s0*X[42] + s1*St_trigDetSumsC::instance()->getPVPDEast();
765 
766  // NoKiller data:
767  // Stored in rich scalers rs12,13,15 for 2011+ data, and available
768  // for the DAQ stream via otherwise empty CTB members of
769  // trigDetSums starting with SL13b
770  X[43] = s0*X[43] + s1*St_trigDetSumsC::instance()->getCTBOrTOFp(); // ZDCXNoKiller
771  X[44] = s0*X[44] + s1*St_trigDetSumsC::instance()->getCTBWest(); // ZDCWestNoKiller
772  X[45] = s0*X[45] + s1*St_trigDetSumsC::instance()->getCTBEast(); // ZDCEastNoKiller
773 
774  // StMagUtilities distortion correction parameters
775  X[46] = s0*X[46] + s1*m_ExB->GetConst_0();
776  X[47] = s0*X[47] + s1*m_ExB->GetConst_1();
777 
778  // SpaceCharge east/west asymmetry
779  X[48] = scE;
780  X[49] = scW;
781  X[50] = s0*X[50] + s1*lastEWRatio*runinfo->spaceCharge();
782 
783  for (i=0;i<24;i++) {
784  X[51+3*i] = scS[i];
785  X[52+3*i] = gapZfitslopeS[i];
786  X[53+3*i] = gapZfitinterceptS[i];
787  }
788 
789  // In calib mode, only fill when doReset (we found an sc)
790  if (doReset || !Calibmode) ntup->Fill(X);
791 
792  if (doReset) {ntent = 0.0; nttrk = 0.0; }
793 
794  }
795 
796  lasttime = thistime;
797 
798  return result;
799 }
800 //_____________________________________________________________________________
801 Int_t StSpaceChargeEbyEMaker::DecideSpaceCharge(int time) {
802 
803  // curhist has only this event
804  // curhist-1 has past two events...
805  // curhist-x has past (x+1) events
806  // curhist-(SCHN-1) == curhist+1 has past SCHN events
807 
808  Bool_t QAout = QAmode || PrePassmode;
809  Bool_t do_auto = kTRUE;
810  Bool_t few_stats = kTRUE;
811  Bool_t large_err = kTRUE;
812  Bool_t large_dif = kTRUE;
813 
814  // Cuts on difference from previous sc measure:
815  // If last event was auto, use MAXDIFFA, else use between MAXDIFFE & MAXDIFFA
816  // scaled by oldness of previous sc measure (curhist-1)
817  float maxdiff,dif=0;
818  if (did_auto)
819  maxdiff = MAXDIFFA;
820  else
821  maxdiff = MAXDIFFA - (MAXDIFFA-MAXDIFFE)*oldness(imodHN(curhist-1));
822 
823  // More than 30 seconds since last used event? Forget it...
824  int timedif = time-lasttime;
825  if (QAout) {
826  gMessMgr->Info() << "time since last event = "
827  << timedif << endm;
828  gMessMgr->Info() << "curhist = "
829  << curhist << endm;
830  }
831  float ntrkstot = 0; // running sum using oldness scale factor
832  float ntrkstotE = 0;
833  float ntrkstotW = 0;
834  Bool_t decideFromData = ((PrePassmode) || (Calibmode) || (lasttime==0) || (timedif < 30));
835  if (decideFromData) {
836 
837  int isc;
838  static int iscMax = 1; // use only one hist for calib mode, and...
839  if (!Calibmode && iscMax<SCHN) iscMax = curhist+1; // don't use uninitialized
840  for (isc=0; isc<iscMax; isc++) {
841  int i = imodHN(curhist-isc);
842  ntrkstot += ntrks[i] * oldness(i);
843  ntrkstotE += ntrksE[i] * oldness(i);
844  ntrkstotW += ntrksW[i] * oldness(i);
845  if (QAout) {
846  if (!isc) gMessMgr->Info("Building with: i, ni, oi, nt:");
847  gMessMgr->Info() << "Building with: " << i << ", "
848  << ntrks[i] << ", " << oldness(i) << ", " << ntrkstot << endm;
849  gMessMgr->Info() << "....east with: " << i << ", "
850  << ntrksE[i] << ", " << oldness(i) << ", " << ntrkstotE << endm;
851  gMessMgr->Info() << "....west with: " << i << ", "
852  << ntrksW[i] << ", " << oldness(i) << ", " << ntrkstotW << endm;
853  }
854 
855  // Too little data collected? Keep trying...
856  few_stats = (IAttr("EastOff") ?
857  (ntrkstotW < MINTRACKS) : // west-only
858  (IAttr("WestOff") ?
859  (ntrkstotE < MINTRACKS) : // east-only
860  (Asymmode ?
861  (ntrkstotE < MINTRACKS || ntrkstotW < MINTRACKS) : // all, Asymmode
862  (ntrkstot < MINTRACKS) ) ) ); // all, not Asymmode
863  if (!few_stats) {
864  BuildHist(i,schist ,schists );
865  BuildHist(i,schistE,schistsE);
866  BuildHist(i,schistW,schistsW);
867  FindSpaceCharge(schist ,sc ,esc );
868  FindSpaceCharge(schistE,scE,escE);
869  FindSpaceCharge(schistW,scW,escW);
870  if (doSecGaps) for (int j=0;j<24;j++) FindSpaceCharge(schistS[j],scS[j],escS[j]);
871  if (QAout) {
872  gMessMgr->Info() << "sc = " << sc << " +/- " << esc << endm;
873  gMessMgr->Info() << "scE = " << scE << " +/- " << escE << endm;
874  gMessMgr->Info() << "scW = " << scW << " +/- " << escW << endm;
875  if (doSecGaps) for (int j=0;j<24;j++) {
876  gMessMgr->Info() << "sc" << Form("%2d",j+1) << " = " << scS[j] << " +/- " << escS[j] << endm;
877  }
878  }
879  large_err = (esc == 0) || (esc > SCALER_ERROR);
880  if (!large_err) {
881  if (PrePassmode) { do_auto=kFALSE; break; }
882  // otherwise, check for big jumps
883  dif = TMath::Abs(sc-lastsc);
884  large_dif = dif > maxdiff;
885  if (!large_dif || Calibmode) {
886  oldevt = evts[i];
887  do_auto=kFALSE;
888  break;
889  }
890  }
891  }
892 
893  // shouldn't need to go past oldest event previously used?
894  // tough to know - allow for now as of 11 Jan 2008
895  // if (evts[i] <= oldevt) break;
896  }
897  if (QAout && (isc == SCHN)) gMessMgr->Info()
898  << "STORED DATA EXHAUSTED: "
899  << SCHN << " events" << endm;
900  }
901 
902  did_auto = do_auto;
903 
904  // In normal mode, do_auto decides whether to use automatic SC from DB
905  // In PrePass mode, do_auto decides when we're ready to stop
906  // In Calib mode, do_auto decides when to save entries and reset
907 
908  if (do_auto) {
909  if (QAout && decideFromData) {
910  if (few_stats) gMessMgr->Info()
911  << "(RECENT) STATS TOO FEW: "
912  << ntrkstot << " / " << ntrkstotE << " / " << ntrkstotW
913  << " (" << MINTRACKS << ")" << endm;
914  else if (large_err) gMessMgr->Info()
915  << "FIT ERROR TOO LARGE: "
916  << esc << " (" << SCALER_ERROR << ")" << endm;
917  else if (large_dif) gMessMgr->Info()
918  << "DIFFERENCE TOO LARGE: "
919  << dif << " (" << maxdiff << ")" << endm;
920  }
921  gMessMgr->Info("using auto SpaceCharge");
922  if (Calibmode) doReset = kFALSE;
923  else m_ExB->AutoSpaceChargeR2();
924  } else {
925  gMessMgr->Info() << "using SpaceCharge = "
926  << sc << " +/- " << esc << " (" << time << ")" << endm;
927  scehist->SetBinContent(evt,sc);
928  scehist->SetBinError(evt,esc);
929  if (PrePassmode) {
930  PrePassdone = kTRUE;
931  return kStStop; // We're happy! Let's stop!
932  }
933  if (Calibmode) {
934  doReset = kTRUE;
935  return kStStop; // We're happy? Let's stop!
936  }
937  else m_ExB->ManualSpaceChargeR2(sc,lastEWRatio);
938  }
939  return kStOk;
940 }
941 //_____________________________________________________________________________
942 void StSpaceChargeEbyEMaker::FindSpaceCharge(TH1F* aschist, float& asc, float& aesc) {
943  aesc = 0.;
944  double res = FindPeak(aschist,aesc);
945  asc = (res > -500. ? res : 0.0);
946 }
947 //_____________________________________________________________________________
948 double StSpaceChargeEbyEMaker::FindPeak(TH1* hist,float& pkwidth) {
949 
950  if (!hist) return -996.;
951  pkwidth = 0.;
952  if (hist->Integral() < 100.0) return -997.;
953  double mn = hist->GetMean();
954  double rms = TMath::Abs(hist->GetRMS());
955  double range = hist->GetXaxis()->GetXmax() - hist->GetXaxis()->GetXmin();
956  mn *= 1.0 + (rms/range); rms *= 1.5;
957  double lr = mn-rms;
958  double ur = mn+rms;
959  double pmax = TMath::Max(hist->GetMaximum(),0.);
960  double lp = pmax*0.001;
961  double up = pmax*10.0;
962  double lw = rms*0.001;
963  double uw = rms*10.0;
964  ga1.SetParameters(pmax,mn,rms*0.5);
965  ga1.SetRange(lr,ur);
966  ga1.SetParLimits(0,lp,up); // Loglikelihood only works with positive functions
967  ga1.SetParLimits(2,lw,uw); // To help the fit
968  hist->Sumw2();
969  int fitResult = hist->Fit(&ga1,
970  (gROOT->GetVersionInt() >= 53000 ? "WLRB0Q" : "LLRB0Q")); // Loglikelihood options changed!
971  hist->Sumw2(kFALSE);
972  ga1.ReleaseParameter(0);
973  ga1.ReleaseParameter(2);
974  if (fitResult != 0) return -999.;
975  double rp = ga1.GetParameter(0);
976  if (rp == lp || rp == up) return -998;
977  pkwidth = TMath::Abs(ga1.GetParError(1));
978  return ga1.GetParameter(1);
979 
980 }
981 //_____________________________________________________________________________
982 void StSpaceChargeEbyEMaker::InitQAHists() {
983 
984  scehist = new TH1F("SpcChgEvt","SpaceCharge fit vs. Event",
985  EVN,0.,EVN);
986  timehist = new TH1F("EvtTime","Event Times",
987  EVN,0.,EVN);
988  dcehist = new TH2F("DcaEve","psDCA vs. Event",
989  EVN,0.,EVN,DCN,DCL,DCH);
990  dcphist = new TH3F("DcaPhi","psDCA vs. Phi",
991  PHN,0,PI2,DCN,DCL,DCH,(QAmode ? 3 : 1),-1.5,1.5);
992 
993  AddHist(scehist);
994  AddHist(timehist);
995  AddHist(dcehist);
996  AddHist(dcphist);
997 
998  if (QAmode) {
999  myhist = new TH3F("SpcEvt","SpaceCharge vs. Phi vs. Event",
1000  EVN,0.,EVN,PHN,0,PI2,SCN2,SCL,SCH);
1001  dcahist = new TH3F("DcaEvt","psDCA vs. Phi vs. Event",
1002  EVN,0.,EVN,PHN,0,PI2,DCN,DCL,DCH);
1003  dczhist = new TH2F("DcaZ","psDCA vs. Z",
1004  //80,-200,200,250,-5.0,5.0);
1005  ZN,ZL,ZH,DCN,DCL,DCH);
1006  myhistN = new TH3F("SpcEvtN","SpaceCharge vs. Phi vs. Event Neg",
1007  EVN,0.,EVN,PHN,0,PI2,SCN2,SCL,SCH);
1008  myhistP = new TH3F("SpcEvtP","SpaceCharge vs. Phi vs. Event Pos",
1009  EVN,0.,EVN,PHN,0,PI2,SCN2,SCL,SCH);
1010  myhistE = new TH3F("SpcEvtE","SpaceCharge vs. Phi vs. Event East",
1011  EVN,0.,EVN,PHN,0,PI2,SCN2,SCL,SCH);
1012  myhistW = new TH3F("SpcEvtW","SpaceCharge vs. Phi vs. Event West",
1013  EVN,0.,EVN,PHN,0,PI2,SCN2,SCL,SCH);
1014  dcahistN = new TH3F("DcaEvtN","psDCA vs. Phi vs. Event Neg",
1015  EVN,0.,EVN,PHN,0,PI2,DCN,DCL,DCH);
1016  dcahistP = new TH3F("DcaEvtP","psDCA vs. Phi vs. Event Pos",
1017  EVN,0.,EVN,PHN,0,PI2,DCN,DCL,DCH);
1018  dcahistE = new TH3F("DcaEvtE","psDCA vs. Phi vs. Event East",
1019  EVN,0.,EVN,PHN,0,PI2,DCN,DCL,DCH);
1020  dcahistW = new TH3F("DcaEvtW","psDCA vs. Phi vs. Event West",
1021  EVN,0.,EVN,PHN,0,PI2,DCN,DCL,DCH);
1022  cutshist = new TH1I("CutsHist","Step at which tracks and events are cut",
1023  64,-0.5,63.5); // Use < 16 for event-wise cuts
1024  AddHist(myhist);
1025  AddHist(dcahist);
1026  AddHist(dczhist);
1027  AddHist(myhistN);
1028  AddHist(myhistP);
1029  AddHist(myhistE);
1030  AddHist(myhistW);
1031  AddHist(dcahistN);
1032  AddHist(dcahistP);
1033  AddHist(dcahistE);
1034  AddHist(dcahistW);
1035  AddHist(cutshist);
1036  }
1037 
1038  if (doNtuple) ntup = new TNtuple("SC","Space Charge",
1039  "sc:dca:zdcx:zdcw:zdce:bbcx:bbcw:bbce:bbcbb:bbcyb:intb:inty:fill:mag:run:event:dcan:dcap:dcae:dcaw:gapf:gapi:gapd:gapfn:gapin:gapdn:gapfp:gapip:gapdp:gapfe:gapie:gapde:gapfw:gapiw:gapdw:usc:uscmode:ugl:zdcc:bbcc:vpdx:vpdw:vpde:zdcxnk:zdcwnk:zdcenk:const0:const1:sce:scw:usce:sc1:gapf1:gapi1:sc2:gapf2:gapi2:sc3:gapf3:gapi3:sc4:gapf4:gapi4:sc5:gapf5:gapi5:sc6:gapf6:gapi6:sc7:gapf7:gapi7:sc8:gapf8:gapi8:sc9:gapf9:gapi9:sc10:gapf10:gapi10:sc11:gapf11:gapi11:sc12:gapf12:gapi12:sc13:gapf13:gapi13:sc14:gapf14:gapi14:sc15:gapf15:gapi15:sc16:gapf16:gapi16:sc17:gapf17:gapi17:sc18:gapf18:gapi18:sc19:gapf19:gapi19:sc20:gapf20:gapi20:sc21:gapf21:gapi21:sc22:gapf22:gapi22:sc23:gapf23:gapi23:sc24:gapf24:gapi24");
1040 
1041  if (doGaps) {
1042  gapZhist = new TH2F("Gaps","Gaps",GZN,GZL,GZH,GN,GL,GH);
1043  gapZhistneg = new TH2F("Gapsneg","Gaps Neg",GZN,GZL,GZH,GN,GL,GH);
1044  gapZhistpos = new TH2F("Gapspos","Gaps Pos",GZN,GZL,GZH,GN,GL,GH);
1045  if (doSecGaps) for (int i=0; i<24; i++) {
1046  int j = i+1;
1047  schistS[i] = new TH1F(Form("SpCh%02d",j),Form("Space Charge%02d",j),SCN1,SCL,SCH);
1048  gapZhistS[i] = new TH2F(Form("Gaps%02d",j),Form("Gaps%02d",j),GZN,GZL,GZH,GN,GL,GH);
1049  gapZhistnegS[i] = new TH2F(Form("Gapsneg%02d",j),Form("Gaps Neg%02d",j),GZN,GZL,GZH,GN,GL,GH);
1050  gapZhistposS[i] = new TH2F(Form("Gapspos%02d",j),Form("Gaps Pos%02d",j),GZN,GZL,GZH,GN,GL,GH);
1051  }
1052  }
1053 
1054 }
1055 //_____________________________________________________________________________
1056 void StSpaceChargeEbyEMaker::WriteQAHists() {
1057 // Only if QAmode or doNtuple or doGaps
1058 
1059  if (!(QAmode || doNtuple || doGaps)) return;
1060 
1061  if (runid == 0) {
1062  gMessMgr->Error("StSpaceChargeEbyEMaker: No runid => No output");
1063  return;
1064  }
1065 
1066  const char* f1 = GetName();
1067  runid -= 1000000*(runid/1000000); // Remove the year, optional
1068  runid *= 100; // limits one run to 100 files!
1069  TString fname = "./hists";
1070  if (PrePassmode) fname.Append("Pre");
1071  gSystem->cd(fname.Data());
1072  while (f1) {
1073  fname = Form("Hist%d.root",runid);
1074  f1 = gSystem->Which(".",fname.Data());
1075  runid++;
1076  delete [] f1;
1077  }
1078 
1079  TFile ff(fname.Data(),"RECREATE");
1080  if (QAmode) {
1081  myhist->Write();
1082  dcehist->Write();
1083  dcphist->Write();
1084  dcahist->Write();
1085  dczhist->Write();
1086  myhistN->Write();
1087  dcahistN->Write();
1088  myhistP->Write();
1089  dcahistP->Write();
1090  myhistE->Write();
1091  dcahistE->Write();
1092  myhistW->Write();
1093  dcahistW->Write();
1094  scehist->Write();
1095  timehist->Write();
1096  cutshist->Write();
1097  }
1098  if (doGaps) {
1099  gapZhist->Write();
1100  gapZhistneg->Write();
1101  gapZhistpos->Write();
1102  if (doSecGaps) for (int i=0; i<24; i++) {
1103  schistS[i]->Write();
1104  gapZhistS[i]->Write();
1105  gapZhistnegS[i]->Write();
1106  gapZhistposS[i]->Write();
1107  }
1108  }
1109  if (doNtuple) ntup->Write();
1110  if (SCcorrection) {
1111  SCcorrection->Write();
1112  SCEWRatio->Write();
1113  GLcorrection->Write();
1114  }
1115  ff.Close();
1116 
1117  gMessMgr->Info() << "QA hists file: " << fname.Data() << endm;
1118 
1119  gSystem->cd("..");
1120 
1121 }
1122 //_____________________________________________________________________________
1123 void StSpaceChargeEbyEMaker::FillQAHists(float DCA, float space, float spaceEW,
1124  int ch, StPhysicalHelixD& hh, int e_or_w) {
1125  // Find a "Phi" for the track
1126  pairD pls = hh.pathLength(97.0);
1127  double pl = pls.first;
1128  if (TMath::Abs(pls.second) < TMath::Abs(pls.first)) pl = pls.second;
1129  StThreeVectorD hh_at_pl = hh.at(pl);
1130  float Phi = hh_at_pl.phi();
1131  while (Phi < 0) Phi += PI2;
1132  while (Phi >= TMath::TwoPi()) Phi -= PI2;
1133 
1134  // To pile all sectors atop each other:
1135  // while (Phi >= TMath::Pi()/6.) Phi -= TMath::Pi()/6.;
1136 
1137  float evtn = ((float) evt) - 0.5;
1138  dcehist->Fill(evtn,DCA);
1139  dcphist->Fill(Phi,DCA,(float) e_or_w);
1140 
1141  if (QAmode) {
1142  myhist->Fill(evtn,Phi,space);
1143  dcahist->Fill(evtn,Phi,DCA);
1144  if (ch > 0) {
1145  myhistP->Fill(evtn,Phi,space);
1146  dcahistP->Fill(evtn,Phi,DCA);
1147  } else {
1148  myhistN->Fill(evtn,Phi,space);
1149  dcahistN->Fill(evtn,Phi,DCA);
1150  }
1151  if (e_or_w > 0) {
1152  myhistW->Fill(evtn,Phi,spaceEW);
1153  dcahistW->Fill(evtn,Phi,DCA);
1154  } else if (e_or_w < 0) {
1155  myhistE->Fill(evtn,Phi,spaceEW);
1156  dcahistE->Fill(evtn,Phi,DCA);
1157  }
1158  if ((e_or_w != 0) && (TMath::Abs(hh.dipAngle()) < 0.05)) dczhist->Fill(hh_at_pl.z(),DCA);
1159  }
1160 }
1161 //_____________________________________________________________________________
1162 int StSpaceChargeEbyEMaker::imodHN(int i) {
1163  // Keep index in bounds of circular queue
1164  return ( i >= SCHN ? imodHN(i-SCHN) : (i < 0 ? imodHN(i+SCHN) : i) );
1165 }
1166 //_____________________________________________________________________________
1167 float StSpaceChargeEbyEMaker::oldness(int i, int j) {
1168  // Deterime how to treat relative "age" of event
1169  // In PrePassmode, earliest events are most important!
1170  float s = 1.0;
1171  if (!PrePassmode) { // Weight newest the most (or evenly for PrePass)
1172  if (j<0) j = curhist;
1173 
1174  // Weight in sub-second intervals by # of events because
1175  // times have only 1 second granularity (best we can do).
1176  // Method assumes time-ordering of events, which is violated
1177  // perhaps only occasionally from DAQ.
1178  int k = i;
1179  while (k!=j) {
1180  k = imodHN(k+1);
1181  if (times[k] != times[i]) { k = imodHN(k-1); break; }
1182  }
1183  // # seconds + fraction of a second:
1184  float time_factor = (times[j]-times[i]) + (1.-(evtstbin[i]/evtstbin[k]));
1185  //float time_factor = (times[j]-times[i]);
1186  float decay_const = -0.12;
1187  // float decay_const = -0.15;
1188  s = exp( decay_const * time_factor );
1189  }
1190  return s;
1191 }
1192 //_____________________________________________________________________________
1193 void StSpaceChargeEbyEMaker::BuildHist(int i, TH1F* aschist, TH1F** aschists) {
1194  // Build up one histogram from several events
1195  aschist->Reset();
1196  int isc = curhist;
1197  aschist->Add(aschists[isc],1.0);
1198  while (isc != i) {
1199  isc = imodHN(isc-1);
1200  aschist->Add(aschists[isc],oldness(isc));
1201  }
1202 }
1203 //_____________________________________________________________________________
1204 void StSpaceChargeEbyEMaker::SetTableName() {
1205  // Problem caused if first event comes later in time than other events.
1206  // Solution: subtract 10 seconds...
1207  TUnixTime ux(GetDateTime(),1);
1208  ux+=-10;
1209  int date,time;
1210  ux.GetGTime(date,time);
1211  gMessMgr->Info() << "first event date = " << date << endm;
1212  gMessMgr->Info() << "first event time = " << time << endm;
1213  tabname = Form("./StarDb/Calibrations/rich/spaceChargeCorR2.%08d.%06d.C",date,time);
1214 
1215  // Set default parameters based on data time
1216  // for non-calib modes
1217  if (TrackInfomode) {
1218  if (TrackInfomode<2) return; // use standard defaults
1219  // ...else show pile-up tracks too
1220  if (TrackInfomode>2) setMinTpcHits(0);
1221  setVtxEmcMatch(0);
1222  setVtxTofMatch(0);
1223  setVtxMinTrks(0);
1224  setReqEmcOrTofMatch(kFALSE);
1225  } else if (date < 20071000) {
1226  setVtxEmcMatch(0);
1227  setVtxTofMatch(0);
1228  setVtxMinTrks(10);
1229  setReqEmcOrTofMatch(kFALSE);
1230  MINTRACKS = 1500;
1231  setVtxVpdAgree(-5);
1232  } else if (date < 20090000) {
1233  setVtxEmcMatch(1);
1234  setVtxTofMatch(0);
1235  setVtxMinTrks(5);
1236  setReqEmcOrTofMatch(kFALSE);
1237  MINTRACKS = 1500;
1238  setVtxVpdAgree(-5);
1239  } else if (date < 20100000) {
1240  setVtxVpdAgree(-5);
1241  }
1242 }
1243 //_____________________________________________________________________________
1244 void StSpaceChargeEbyEMaker::WriteTableToFile(){
1245  gMessMgr->Info() << "Writing new table to:\n "
1246  << tabname.Data() << endm;
1247  TString dirname = gSystem->DirName(tabname.Data());
1248  TString estr = dirname;
1249  estr.Prepend("mkdir -p ");
1250  gSystem->Exec(estr.Data());
1251  if (gSystem->OpenDirectory(dirname.Data())==0) {
1252  if (gSystem->mkdir(dirname.Data())) {
1253  gMessMgr->Warning() << "Directory creation failed for:\n " << dirname
1254  << "\n Putting table file in current directory" << endm;
1255  tabname.Remove(0,tabname.Last('/')).Prepend(".");
1256  }
1257  }
1258  ofstream *out = new ofstream(tabname.Data());
1259  SCTable()->SavePrimitive(*out,"");
1260  return;
1261 }
1262 //_____________________________________________________________________________
1263 St_spaceChargeCor* StSpaceChargeEbyEMaker::SCTable() {
1264  St_spaceChargeCor* table = new St_spaceChargeCor("spaceChargeCorR2",1);
1265  spaceChargeCor_st* row = table->GetTable();
1266  row->fullFieldB = sc;
1267  row->halfFieldB = sc;
1268  row->zeroField = (float) evt;
1269  row->halfFieldA = sc;
1270  row->fullFieldA = sc;
1271  row->satRate = 1.0;
1272  row->factor = 1.0;
1273  row->detector = 3;
1274  row->offset = 0;
1275  row->ewratio = m_ExB->CurrentSpaceChargeEWRatio();
1276  table->SetNRows(1);
1277  return table;
1278 }
1279 //_____________________________________________________________________________
1280 float StSpaceChargeEbyEMaker::FakeAutoSpaceCharge() {
1281  // Use this to mimic what is done with the scalers, using your own code
1282  float zdcsum = runinfo->zdcWestRate()+runinfo->zdcEastRate();
1283  float sc = 6e-8 * zdcsum;
1284 
1285  return sc;
1286 }
1287 //_____________________________________________________________________________
1288 void StSpaceChargeEbyEMaker::FillGapHists(StTrack* tri, StPhysicalHelixD& hh,
1289  int e_or_w, int ch) {
1290  St_tpcPadConfigC* pads = St_tpcPadConfigC::instance();
1291  float fsign = ( event->runInfo()->magneticField() < 0 ? -1 : 1 );
1292  StPtrVecHit hts = tri->detectorInfo()->hits(kTpcId);
1293  float gap = 0.; float zgap = 0.; int ct=0;
1294  float GAPRADIUS = 121.8; // messy to get from the database (see StMagUtilities)
1295  Int_t sector = -1;
1296  for (UInt_t ht=0; ht<hts.size(); ht++) {
1297  StTpcHit* hit = (StTpcHit*) hts[ht];
1298  Int_t prow = hit->padrow();
1299  Int_t sec = hit->sector();
1300  int lastInner = pads->innerPadRows(sec);
1301  if ((prow != lastInner) && (prow != lastInner+1)) continue;
1302  float gsign = ( prow > lastInner ? -1 : 1 );
1303  const StThreeVectorF& hp = hit->position();
1304 
1305  // Avoid sector edges
1306  float hphi = hp.phi() + TMath::TwoPi();
1307  while (hphi > TMath::Pi()/12.) hphi -= TMath::Pi()/6.;
1308  if (TMath::Abs(hphi) > 0.75*TMath::Pi()/12.) break;
1309 
1310  float distToOuterRow = pads->radialDistanceAtRow(sec,lastInner+1) - GAPRADIUS;
1311  float distToInnerRow = GAPRADIUS - pads->radialDistanceAtRow(sec,lastInner + inGapRow - 13);
1312  zgap += (hp.z() / (distToInnerRow+distToOuterRow)) *
1313  ( prow == lastInner ? distToOuterRow : distToInnerRow ); // ~z at gap
1314 
1315  // Measurement method described at:
1316  // http://drupal.star.bnl.gov/STAR/blog/genevb/2010/feb/21/gridleak-update-using-residuals-along-padrows
1317  Double_t residual = hh.geometricSignedDistance(hp.x(),hp.y());
1318 
1319  sector = ( sec == sector || sector < 0 ? sec : 0 ); // 0 if crossing sectors
1320  if (sector == 0) return; // don't use sector-crossing tracks!
1321  Double_t sector_angle = (TMath::Pi()/6.) * (sec < 13 ? 3 - sec : sec - 21);
1322  Double_t pathlen = hh.pathLength(hp.x(),hp.y());
1323  Double_t theta = TMath::ATan2(hh.cy(pathlen),hh.cx(pathlen)) - sector_angle;
1324  Double_t phi = hp.phi() - sector_angle;
1325 
1326  Double_t Eff1 = TMath::Cos(theta);
1327  Double_t Eff2 = 1;
1328  Double_t Eff3 = 0;
1329 
1330  Float_t x1[3],x2[3];
1331  x1[0] = hp.perp()*TMath::Sin(phi);
1332  x1[1] = hp.perp()*TMath::Cos(phi);
1333  x1[2] = hp.z();
1334  m_ExB->UndoGridLeakDistortion(x1,x2,sec);
1335  Double_t dX = x2[0]-x1[0];
1336  if (TMath::Abs(dX) > 1e-20) {
1337  // warning: no available GridLeak calculation may lead
1338  // to slightly different results
1339  Eff3 = gsign * ((x2[1]-x1[1])/dX) * TMath::Sin(theta);
1340  x1[0] = 0;
1341  m_ExB->UndoGridLeakDistortion(x1,x2,sec);
1342  Eff2 = (x2[0]-x1[0])/dX;
1343  }
1344 
1345  Double_t DistortionX = Eff2 * residual / (Eff1 + Eff3);
1346 
1347  gap += fsign * gsign * DistortionX;
1348  ct++;
1349  }
1350 
1351  sector = (doSecGaps ? sector - 1 : -1);
1352  float abs_zgap = TMath::Abs(zgap);
1353  if ((ct==2) && (abs_zgap<200.0) && (abs_zgap>10.0)) {
1354  gapZhist->Fill(zgap,gap);
1355  if (ch==1) gapZhistpos->Fill(zgap,gap);
1356  else gapZhistneg->Fill(zgap,gap);
1357  if (sector >= 0) {
1358  gapZhistS[sector]->Fill(zgap,gap);
1359  if (ch==1) gapZhistposS[sector]->Fill(zgap,gap);
1360  else gapZhistnegS[sector]->Fill(zgap,gap);
1361  }
1362 
1363  if (abs_zgap<150 && abs_zgap>25) { // Restrict the Z range further
1364  // normalize at z=100cm from ggrid
1365  float gap_scaled = (gap * 100.0) / (ZGGRID - abs_zgap);
1366  //float gap_scaled = (gap * 100.0) / (350.0 - abs_zgap);
1367  float z_beyond = ZGGRID+1.0;
1368  gapZhist->Fill(z_beyond,gap_scaled);
1369  if (ch==1) gapZhistpos->Fill(z_beyond,gap_scaled);
1370  else gapZhistneg->Fill(z_beyond,gap_scaled);
1371  if (sector >= 0) {
1372  gapZhistS[sector]->Fill(z_beyond,gap_scaled);
1373  if (ch==1) gapZhistposS[sector]->Fill(z_beyond,gap_scaled);
1374  else gapZhistnegS[sector]->Fill(z_beyond,gap_scaled);
1375  }
1376  }
1377  }
1378 }
1379 //_____________________________________________________________________________
1380 void StSpaceChargeEbyEMaker::DetermineGaps() {
1381  TString GapStr = "Gap measurements: fit slope & intercept, div slope";
1382 
1383  if (doSecGaps) {
1384  GapStr += "\nneg (by sector):";
1385  for (int i=0; i<24; i++) {
1386  GapStr += Form("\n%2d : ",i+1);
1387  GapStr += DetermineGapHelper(gapZhistnegS[i],gapZfitslopenegS[i],gapZfitinterceptnegS[i],gapZdivslopenegS[i]);
1388  }
1389  GapStr += "\npos (by sector):";
1390  for (int i=0; i<24; i++) {
1391  GapStr += Form("\n%2d : ",i+1);
1392  GapStr += DetermineGapHelper(gapZhistposS[i],gapZfitslopeposS[i],gapZfitinterceptposS[i],gapZdivslopeposS[i]);
1393  }
1394  GapStr += "\nall (by sector):";
1395  for (int i=0; i<24; i++) {
1396  GapStr += Form("\n%2d : ",i+1);
1397  GapStr += DetermineGapHelper(gapZhistS[i],gapZfitslopeS[i],gapZfitinterceptS[i],gapZdivslopeS[i]);
1398  }
1399  GapStr += "\nAll sectors:";
1400  }
1401 
1402  GapStr += "\nneg: ";
1403  GapStr += DetermineGapHelper(gapZhistneg,gapZfitslopeneg,gapZfitinterceptneg,gapZdivslopeneg);
1404  GapStr += "\npos: ";
1405  GapStr += DetermineGapHelper(gapZhistpos,gapZfitslopepos,gapZfitinterceptpos,gapZdivslopepos);
1406  GapStr += "\nall: ";
1407  GapStr += DetermineGapHelper(gapZhist,gapZfitslope,gapZfitintercept,gapZdivslope);
1408 
1409  LOG_INFO << GapStr << endm;
1410 }
1411 //_____________________________________________________________________________
1412 TString StSpaceChargeEbyEMaker::DetermineGapHelper(TH2F* gh,
1413  float& fitslope, float& fitintercept, float& divslope) {
1414 
1415  TString result = "n/a";
1416  if (gh->GetEntries() < 10) { // don't waste time fitting ~empty histograms
1417  fitslope = 0; fitintercept = 0; divslope = 0;
1418  egapZfitslope = 0; egapZfitintercept = 0; egapZdivslope = 0;
1419  return result;
1420  }
1421  ga1.SetParameters(gh->GetEntries()/(16.*2.*10.),0.,0.1);
1422  ga1.SetParLimits(0,0.001,10.0*gh->GetEntries()); // Loglikelihood only works with positive functions
1423  gh->FitSlicesY(&ga1,1,0,20,"LB0Q"); // gapZhist bin contents are integers
1424  ga1.ReleaseParameter(0);
1425  const char* hn = gh->GetName();
1426  TH1D* GapsChi2 = (TH1D*) gDirectory->Get(Form("%s_chi2",hn));
1427  TH1D* GapsAmp = (TH1D*) gDirectory->Get(Form("%s_0",hn));
1428  TH1D* GapsMean = (TH1D*) gDirectory->Get(Form("%s_1",hn));
1429  TH1D* GapsRMS = (TH1D*) gDirectory->Get(Form("%s_2",hn));
1430 
1431  divslope = GapsMean->GetBinContent(GZN);
1432  egapZdivslope = TMath::Abs(GapsMean->GetBinError(GZN));
1433 
1434  // Fit only from |Z| = 25-150cm
1435  // Zero bins -25 thru 25cm (depends upon GZN,GZL,GZH definitions)
1436  GapsMean->SetBinContent(8,0); GapsMean->SetBinError(8,0);
1437  GapsMean->SetBinContent(9,0); GapsMean->SetBinError(9,0);
1438 
1439  ln1.SetRange(-150.,150.);
1440  GapsMean->Fit(&ln1,"R0Q");
1441  fitslope = ln1.GetParameter(1);
1442  fitintercept = ln1.GetParameter(0);
1443 
1444  // remaining variables save the last histogram fit: "Gaps" for now
1445 
1446  egapZfitslope = TMath::Abs(ln1.GetParError(1));
1447  egapZfitintercept = TMath::Abs(ln1.GetParError(0));
1448 
1449  ln1.SetRange(-150.,0.);
1450  GapsMean->Fit(&ln1,"R0Q");
1451  gapZfitslopeeast = ln1.GetParameter(1);
1452  gapZfitintercepteast = ln1.GetParameter(0);
1453 
1454  ln1.SetRange(0.,150.);
1455  GapsMean->Fit(&ln1,"R0Q");
1456  gapZfitslopewest = ln1.GetParameter(1);
1457  gapZfitinterceptwest = ln1.GetParameter(0);
1458 
1459  delete GapsChi2;
1460  delete GapsAmp;
1461  delete GapsMean;
1462  delete GapsRMS;
1463 
1464  result = Form("%7.4f+/-%7.4f , %7.4f+/-%7.4f , %7.4f+/-%7.4f",
1465  fitslope,egapZfitslope,fitintercept,egapZfitintercept,
1466  divslope,egapZdivslope);
1467  return result;
1468 }
1469 //_____________________________________________________________________________
1470 float StSpaceChargeEbyEMaker::EvalCalib(TDirectory* hdir) {
1471 
1472  if (hdir) {
1473  dcehist = static_cast<TH2F*>(hdir->Get("DcaEve"));
1474  timehist = static_cast<TH1F*>(hdir->Get("EvtTime"));
1475  scehist = static_cast<TH1F*>(hdir->Get("SpcChgEvt"));
1476  if (!(dcehist && timehist && scehist)) {
1477  LOG_ERROR << "Problems finding SC histograms" << endm;
1478  return 999.;
1479  }
1480  }
1481 
1482  // Other counts
1483  float spc = (float) (scehist->GetEntries());
1484  float dcc = (float) (dcehist->GetEntries());
1485  float evc = (float) (timehist->GetEntries());
1486 
1487  float hm=0,hw=0,gm=0,gw=0,gm1=0,gw1=0,gme=0,gwe=0,pm=0,pw=0,epsec=0,frac=0,wid=9.99;
1488  TF1* pl0 = 0;
1489 
1490  if (dcc>0) {
1491  TH1D* dcaproj = dcehist->ProjectionY();
1492 
1493  // Initial fits to DCA distribution
1494  ga1.SetParameters(1.,0.,1.);
1495  dcaproj->Fit(&ga1,"Q");
1496  hm = ga1.GetParameter(1);
1497  hw = TMath::Abs(ga1.GetParameter(2));
1498  float hd = 0.6*hw;
1499  ga1.SetRange(hm-hd,hm+hd);
1500  dcaproj->Fit(&ga1,"RQ");
1501  gm = ga1.GetParameter(1);
1502  gw = TMath::Abs(ga1.GetParameter(2));
1503  gm1 = gm;
1504  gw1 = gw;
1505 
1506  // Iterate fit to get best answer
1507  ga1.SetRange(gm-0.9*gw,gm+0.9*gw);
1508  dcaproj->Fit(&ga1,"RQ");
1509  gm = ga1.GetParameter(1);
1510  gw = TMath::Abs(ga1.GetParameter(2));
1511  ga1.SetRange(gm-0.8*gw,gm+0.8*gw);
1512  dcaproj->Fit(&ga1,"RQ");
1513  gm = ga1.GetParameter(1);
1514  gw = TMath::Abs(ga1.GetParameter(2));
1515  ga1.SetRange(gm-0.7*gw,gm+0.7*gw);
1516  dcaproj->Fit(&ga1,"RQ");
1517  gm = ga1.GetParameter(1);
1518  gw = TMath::Abs(ga1.GetParameter(2));
1519  gme = TMath::Abs(ga1.GetParError(1));
1520  gwe = TMath::Abs(ga1.GetParError(2));
1521 
1522  // Quality measures
1523  wid = gw1*gw1+gw*gw;
1524  if (wid>0) wid = TMath::Min(10.,TMath::Log10(wid));
1525  }
1526 
1527  if (spc>0) {
1528  // Average SC
1529  scehist->Fit("pol0","Q");
1530  pl0 = scehist->GetFunction("pol0");
1531  if (pl0) {
1532  pm = pl0->GetParameter(0);
1533  pw = pl0->GetParError(0);
1534  }
1535  }
1536 
1537  if (evc) {
1538  timehist->GetXaxis()->SetRange(1,(int) evc);
1539  timehist->Fit("pol0","LQ");
1540  pl0 = timehist->GetFunction("pol0");
1541  if (pl0) epsec = pl0->GetParameter(0); // events per second
1542 
1543  // Quality measures
1544  frac = spc/evc; // fraction of events for which SC was found
1545  }
1546 
1547 
1548  float code=0;
1549  if (frac<0.2) code = 1. + frac; // code = 1.x
1550  else if (wid>0) code = 2. + 0.1*wid; // code = 2.x
1551 
1552  LOG_INFO << Form("SCeval: %f %f %f %f %f %f %f %f %f %f %f %f %f %f\n",
1553  hm,hw,gm,gw,pm,pw,gm1,gw1,spc,dcc,evc,epsec,gme,gwe) << endm;
1554 
1555  if (code>0) {
1556  LOG_ERROR << "CheckFail: Break of SpaceCharge performance! code = " << code << endm;
1557  }
1558 
1559  return code;
1560 }
1561 //_____________________________________________________________________________
1562 // $Id: StSpaceChargeEbyEMaker.cxx,v 1.72 2018/10/19 21:06:17 genevb Exp $
1563 // $Log: StSpaceChargeEbyEMaker.cxx,v $
1564 // Revision 1.72 2018/10/19 21:06:17 genevb
1565 // Clean up after Y. Fisyak modifications (which were for iTPC, not dE/dx), and use new PredictSpaceCharge() using real hit radii
1566 //
1567 // Revision 1.71 2018/10/17 20:45:27 fisyak
1568 // Restore update for Run XVIII dE/dx calibration removed by Gene on 08/07/2018
1569 //
1570 // Revision 1.70 2018/09/21 18:22:47 genevb
1571 // Bug fix for wrong vertex daughter PCTs function
1572 //
1573 // Revision 1.69 2018/06/08 16:39:59 genevb
1574 // Needs explicity include of TROOT.h
1575 //
1576 // Revision 1.68 2018/06/08 15:35:06 genevb
1577 // Needs explicity include of TSystem.h
1578 //
1579 // Revision 1.67 2018/02/15 03:25:29 genevb
1580 // Restore prepass settings for old data
1581 //
1582 // Revision 1.66 2017/04/12 19:47:25 genevb
1583 // Use generic GridLeak function
1584 //
1585 // Revision 1.65 2017/02/14 23:38:38 fisyak
1586 // Adjustment to structure changes in StTpcdEdxCorrection
1587 //
1588 // Revision 1.64 2015/06/30 21:44:31 genevb
1589 // Use an initialization call for StMagUtilities for each event
1590 //
1591 // Revision 1.63 2015/05/23 04:26:07 genevb
1592 // More vertex selection criteria: PCT daughters, and VPD z agreement
1593 //
1594 // Revision 1.62 2015/05/19 19:36:09 genevb
1595 // Code cleanup in preparation for C++11
1596 //
1597 // Revision 1.61 2015/05/15 14:34:45 genevb
1598 // Fix incorrect memset usage (RT 3093)
1599 //
1600 // Revision 1.60 2015/03/11 21:38:52 genevb
1601 // HFT era: no tracks with PXL or IST hits in calibration
1602 //
1603 // Revision 1.59 2014/11/17 20:49:09 genevb
1604 // Store east and west gapf in the ntuple
1605 //
1606 // Revision 1.58 2014/10/23 21:07:23 genevb
1607 // Add GridLeak-by-sector codes, East/WestOff handling, and some code reformatting
1608 //
1609 // Revision 1.57 2014/07/23 17:58:46 genevb
1610 // Machinery for sector-by-sector Gaps (GridLeak) measurements
1611 //
1612 // Revision 1.56 2014/06/26 22:06:26 fisyak
1613 // New Tpc Alignment, v632
1614 //
1615 // Revision 1.55 2014/05/15 17:14:03 genevb
1616 // Minor tweaks for TrackInfo mode
1617 //
1618 // Revision 1.54 2014/05/02 02:38:07 genevb
1619 // TrackInfo mode with pile-up tracks too
1620 //
1621 // Revision 1.53 2014/01/02 20:54:28 genevb
1622 // TrackInfomode, and Basic E/W asymmetry functionality
1623 //
1624 // Revision 1.52 2013/09/25 20:55:51 genevb
1625 // Allow use of multiple PPVF vertices, introduce EmcOrTofMatch, keep track of Predict...() cuts
1626 //
1627 // Revision 1.51 2013/04/26 20:00:54 genevb
1628 // Protection against 0 entry histos for EvalCalib()
1629 //
1630 // Revision 1.50 2013/03/11 20:04:31 genevb
1631 // make ugl and average over data
1632 //
1633 // Revision 1.49 2013/03/09 23:37:35 genevb
1634 // Add NoKiller ZDC data to ntuple
1635 //
1636 // Revision 1.48 2013/03/07 23:12:30 genevb
1637 // Improve FindPeak(), particularly for sc hists at high lumi, and add StMagUtil:const0,1 to ntuple
1638 //
1639 // Revision 1.47 2013/02/19 21:07:58 genevb
1640 // Print out the sc and GL correction formulas
1641 //
1642 // Revision 1.46 2012/12/28 22:04:25 genevb
1643 // Improve chances of fits succeeding
1644 //
1645 // Revision 1.45 2012/12/15 03:13:50 genevb
1646 // Store used calibrations in histogram files
1647 //
1648 // Revision 1.44 2012/10/01 17:50:07 genevb
1649 // Reduce some overhead DB queries by being more specific about needed tables
1650 //
1651 // Revision 1.43 2012/09/13 20:58:56 fisyak
1652 // Corrections for iTpx
1653 //
1654 // Revision 1.42 2012/08/18 02:11:59 genevb
1655 // Expand SC hist ranges, add VPD to ntuple
1656 //
1657 // Revision 1.41 2012/04/25 19:23:55 genevb
1658 // Pointing angle near vertex needed for improved PredictSpaceCharge
1659 //
1660 // Revision 1.40 2012/01/14 00:21:04 genevb
1661 // Add code for EMC match, set default to required
1662 //
1663 // Revision 1.39 2011/10/27 23:11:03 genevb
1664 // Account for pointing error in sDCA and predicted SpaceCharge
1665 //
1666 // Revision 1.38 2011/10/26 15:33:49 genevb
1667 // Avoid floating point exception in Loglikelihood fits of gaussian peaks
1668 //
1669 // Revision 1.37 2011/04/18 17:36:52 genevb
1670 // Expanded range for sc hists
1671 //
1672 // Revision 1.36 2011/02/10 18:31:45 genevb
1673 // Restore corrected coincidence rates, add QA histogram of where events/tracks are cut
1674 //
1675 // Revision 1.35 2011/02/09 21:56:50 genevb
1676 // Version which can work in SL10k
1677 //
1678 // Revision 1.34 2011/02/09 21:11:36 genevb
1679 // Parameters need to be available in normal event-by-event mode too
1680 //
1681 // Revision 1.33 2011/02/09 16:24:18 genevb
1682 // Allow for historical operating parameters in Prepass mode
1683 //
1684 // Revision 1.32 2010/11/17 17:23:33 genevb
1685 // Include corrected coincidence rates in ntuple
1686 //
1687 // Revision 1.31 2010/07/09 19:01:54 genevb
1688 // Add TOF matching
1689 //
1690 // Revision 1.30 2010/06/09 20:24:53 genevb
1691 // Modify interface to allow EMC and TOF matching requirements (needs implementation)
1692 //
1693 // Revision 1.29 2010/02/25 21:50:18 genevb
1694 // Pass sector number to StMagUtilities, correct unsigned int usage
1695 //
1696 // Revision 1.28 2010/02/24 22:58:21 genevb
1697 // Even the highest ranked vertex might need to be above a minimum ranking
1698 //
1699 // Revision 1.27 2010/02/23 23:59:41 genevb
1700 // Reduce positional biases in GridLeak measurements
1701 //
1702 // Revision 1.26 2010/01/27 15:11:00 fisyak
1703 // eliminate access to StTpcDbMaker, use directly gStTpcDb
1704 //
1705 // Revision 1.25 2010/01/27 14:42:33 fisyak
1706 // Use StMagUtilities::Instance instead of asking tpcHitMoverMaker
1707 //
1708 // Revision 1.24 2009/11/16 22:02:19 genevb
1709 // Loosen nDaughters cut, add BEMCmatch cut, PCT hits cut, enable padrow 13 for Run 9+
1710 //
1711 // Revision 1.23 2009/11/10 20:54:13 fisyak
1712 // pams Cleanup
1713 //
1714 // Revision 1.22 2008/07/24 19:17:55 genevb
1715 // SpaceChargeEWRatio must be written to prepass output table
1716 //
1717 // Revision 1.21 2008/07/21 21:17:10 genevb
1718 // No call to EvalCalib() if Prepass runs successfully
1719 //
1720 // Revision 1.20 2008/07/17 23:35:33 jeromel
1721 // Small format change
1722 //
1723 // Revision 1.19 2008/07/16 03:42:15 genevb
1724 // Use CheckFail in failed performance check error message
1725 //
1726 // Revision 1.18 2008/07/15 22:30:38 genevb
1727 // Added evaluation of calibration performance
1728 //
1729 // Revision 1.17 2008/04/30 14:52:15 genevb
1730 // Reduce pileup contributions
1731 //
1732 // Revision 1.16 2008/01/14 19:22:49 genevb
1733 // Fine tuning of parameters, removal of excess text in messages
1734 //
1735 // Revision 1.15 2007/01/25 19:04:04 perev
1736 // GMT fix
1737 //
1738 // Revision 1.14 2007/01/24 21:42:22 perev
1739 // GMT conversion fixed
1740 //
1741 // Revision 1.13 2006/12/16 01:00:58 genevb
1742 // Better handling of zero magnetic field
1743 //
1744 // Revision 1.12 2006/08/15 23:40:59 genevb
1745 // Averaging was done improperly in DontReset mode
1746 //
1747 // Revision 1.11 2006/07/17 20:13:08 genevb
1748 // Disallow SVT points on tracks
1749 //
1750 // Revision 1.10 2006/07/02 23:22:36 genevb
1751 // Allow for SVT/SSD hits on tracks (necessary for ITTF)
1752 //
1753 // Revision 1.9 2006/06/01 17:27:11 genevb
1754 // Bug fix: gapd and gapf backwards; Improvements: gap fit intercepts, hist and fit ranges
1755 //
1756 // Revision 1.8 2006/01/05 19:12:53 genevb
1757 // Added calib mode
1758 //
1759 // Revision 1.7 2005/12/07 19:45:46 perev
1760 // Histos diconnected from directory. EndCrashFix
1761 //
1762 // Revision 1.6 2005/04/21 19:38:20 genevb
1763 // Additional code for studying SpaceCharge
1764 //
1765 // Revision 1.5 2005/02/16 17:18:06 genevb
1766 // Fill StEvent info on SpaceCharge
1767 //
1768 // Revision 1.4 2004/08/13 20:49:12 genevb
1769 // Improve upon keeping method locked on for each event, and timestamp change
1770 //
1771 // Revision 1.3 2004/08/02 01:19:27 genevb
1772 // minor fixes for getting directories correct
1773 //
1774 // Revision 1.2 2004/07/01 01:46:04 genevb
1775 // Slightly larger margin for full vs half field differences
1776 //
1777 // Revision 1.1 2004/06/30 23:16:00 genevb
1778 // Introduction of StSpaceChargeEbyEMaker
1779 //
1780 //
Bool_t trackOnEmc(StThreeVectorD *position, StThreeVectorD *momentum, const StTrack *const track, Double_t magField, Double_t emcRadius=225.405) const
Track projection utility.
Definition: StHit.h:125
double distance(const StThreeVector< double > &p, bool scanPeriods=true) const
minimal distance between point and helix
Definition: StHelix.cc:240
double cx(double s) const
pointing vector of helix at point s
Definition: StHelix.hh:203
pair< double, double > pathLength(double r) const
path length at given r (cylindrical r)
Definition: StHelix.cc:351
unsigned char matchFlag() const
Matching information.
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
Definition: Stypes.h:42
Definition: Stypes.h:51
Int_t getBin(const Float_t phi, const Float_t eta, Int_t &m, Int_t &e, Int_t &s) const
Definition: StEmcGeom.h:321
Definition: Stypes.h:41
double Dca(const double point[3], double *dcaErr=0) const
DCA to given space point (with error matrix)