Back to home page

sPhenix code displayed by LXR

 
 

    


Warning, file /analysis/LuminosityCounterGoodRuns/src/triggercountmodule.cc was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 #include "triggercountmodule.h"
0002 #include <globalvertex/GlobalVertexMap.h>
0003 #include <globalvertex/GlobalVertexMapv1.h>
0004 #include <globalvertex/GlobalVertex.h>
0005 #include <fun4all/Fun4AllReturnCodes.h>
0006 #include <phool/PHCompositeNode.h>
0007 #include <phool/getClass.h>
0008 #include <ffaobjects/EventHeaderv1.h>
0009 #include <ffarawobjects/Gl1Packetv2.h>
0010 using namespace std;
0011 
0012 
0013 //____________________________________________________________________________..
0014 triggercountmodule::triggercountmodule(const std::string &filename, int rn, int segn, int maxseg, int debug, const std::string &name, int clt):
0015   SubsysReco(name)//).c_str())
0016 {
0017   _rn = rn;
0018   _filename = filename;
0019   _name = name;
0020 
0021   _tree = new TTree("_tree","a tree to count events");
0022   _debug = debug;
0023   _nseg = segn;
0024   _lastseg = maxseg;
0025   _clt = clt;
0026 }
0027 
0028 //____________________________________________________________________________..
0029 triggercountmodule::~triggercountmodule()
0030 {
0031 
0032 }
0033 
0034 //____________________________________________________________________________..
0035 int triggercountmodule::Init(PHCompositeNode *topNode)
0036 {
0037   
0038   if(_debug > 1) cout << "Begin init: " << endl;
0039   
0040   _tree->Branch("badFlag",&_badFlag,"badFlag/I");
0041   _tree->Branch("startLive",_startLive,"startLive[64]/l");
0042   _tree->Branch("endLive",_endLive,"endLive[64]/l");
0043   _tree->Branch("startScal",_startScal,"startScal[64]/l");
0044   _tree->Branch("endScal",_endScal,"endScal[64]/l");
0045   _tree->Branch("nevt",&_evtn,"nevt/l");
0046   _tree->Branch("avgPS",_avgPS,"avgPS[64]/D");
0047   _tree->Branch("trigCounts",_trigCounts,"trigCounts[6][64]/l");
0048   _tree->Branch("eMBDlive",_eMBDlive,"eMBDlive[5]/D");
0049   _tree->Branch("startBCO",&_startBCO,"startBCO/l");
0050   _tree->Branch("endBCO",&_endBCO,"endBCO/l");
0051   _tree->Branch("nBunch",&_nBunch,"nBunch/I");
0052   _tree->Branch("startRaw",_startRaw,"startRaw[64]/l");
0053   _tree->Branch("endRaw",_endRaw,"endRaw[64]/l");
0054   _mbzhist = new TH1D("mbzhist","",600,-300,300);
0055   _bunchHist = new TH1D("bunchHist","",121,0,121);
0056   if(_debug > 0)
0057     {
0058       cout << "Printing tree: " << endl;
0059       _tree->Print();
0060     }
0061   return Fun4AllReturnCodes::EVENT_OK;
0062 }
0063 
0064 //____________________________________________________________________________..
0065 int triggercountmodule::InitRun(PHCompositeNode *topNode)
0066 {
0067   if(_debug > 1) cout << "Initializing!" << endl;
0068   return Fun4AllReturnCodes::EVENT_OK;
0069 }
0070 
0071 //____________________________________________________________________________..
0072 
0073 
0074 int triggercountmodule::process_event(PHCompositeNode *topNode)
0075 {
0076 
0077   if(_debug > 1) cout << endl << endl << endl << "triggercountmodule: Beginning event processing" << endl;
0078 
0079     Gl1Packetv2* gl1 = findNode::getClass<Gl1Packetv2>(topNode, "14001");
0080   if(!gl1)
0081     {
0082       cout << "No trigger info!" << endl;
0083       return Fun4AllReturnCodes::ABORTRUN;
0084     }
0085 
0086   _bunchHist->Fill(gl1->lValue(0,"BunchNumber"));
0087   
0088   if(_evtn == 0)
0089     {
0090       _startBCO = gl1->lValue(0,"BCO");
0091       for(int i=0; i<64; ++i)
0092     {
0093       _startLive[i] = gl1->lValue(i,1);
0094       if(_debug > 2) cout << gl1->lValue(i,1) << endl;
0095       _startScal[i] = gl1->lValue(i,2);
0096       _startRaw[i] = gl1->lValue(i,0);
0097       _prevcheck[0] = gl1->lValue(18,0);
0098       _prevcheck[1] = gl1->lValue(22,0);
0099       _prevcheck[2] = gl1->lValue(18,1);
0100       _prevcheck[3] = gl1->lValue(22,1);
0101     }
0102     }
0103   //cout << _evtn << endl;
0104 
0105   _endBCO = gl1->lValue(0,"BCO");
0106   for(int i=0; i<64; ++i)
0107     {
0108       _endLive[i] = gl1->lValue(i,1);
0109       if(_debug > 2) cout << "end" << gl1->lValue(i,1) << endl;
0110       _endScal[i] = gl1->lValue(i,2);
0111       _endRaw[i] = gl1->lValue(i,0);
0112     }
0113 
0114   float ld18 = _endLive[18] - _prevcheck[2];
0115   float ld22 = _endLive[22] - _prevcheck[3];
0116   float rd18 = _endRaw[18] - _prevcheck[0];
0117   float rd22 = _endRaw[22] - _prevcheck[1];
0118 
0119   float lt18 = ld18/rd18;
0120   float lt22 = ld22/rd22;
0121 
0122   _isblt = 0;
0123   if(_evtn != 0 && _clt)
0124     {
0125       if(lt18 < 0.1 || lt22 < 0.1)
0126     {
0127       _isblt = 1;
0128     }
0129     }
0130   
0131   if(_debug > 1) cout << "Getting gl1 trigger vector from: " << gl1 << endl;
0132   ++_evtn;
0133   long long unsigned int _trigvec = gl1->getScaledVector();
0134 
0135   int ismbtrig = (_trigvec >> 10) & 1;
0136 
0137   GlobalVertexMap* gvtxmap = findNode::getClass<GlobalVertexMapv1>(topNode, "GlobalVertexMap");
0138 
0139   float zvtx = NAN;
0140 
0141   if(gvtxmap)
0142     {
0143       GlobalVertex* gvtx = gvtxmap->begin()->second;
0144       if(gvtx)
0145         {
0146           auto startIter = gvtx->find_vertexes(_vtxtype);
0147           auto endIter = gvtx->end_vertexes();
0148           for(auto iter = startIter; iter != endIter; ++iter)
0149             {
0150               const auto &[type, vertexVec] = *iter;
0151               if(type != _vtxtype) continue;
0152               for(const auto *vertex : vertexVec)
0153                 {
0154                   if(!vertex) continue;
0155                   zvtx = vertex->get_z();
0156                 }
0157             }
0158         }
0159     }
0160   
0161   int zthresh[5] = {30,60,200,10,1000};
0162 
0163 
0164 
0165   for(int i=0; i<5; ++i)
0166     {
0167       if(abs(zvtx) < zthresh[i] && !std::isnan(zvtx))
0168     {
0169       for(int j=0; j<64; ++j)
0170         {
0171           if((_trigvec >> j) & 1)
0172         {
0173           if(!_isblt) _trigCounts[i][j]++;
0174         }
0175         }
0176     }
0177     }
0178 
0179   for(int j=0; j<64; ++j)
0180     {
0181       if(((_trigvec >> j) & 1) && !_isblt) _trigCounts[5][j]++;
0182     }
0183 
0184   if(!ismbtrig)
0185     {
0186       if(_debug > 2) cout << "non MB event, not counting" << endl;
0187       return Fun4AllReturnCodes::ABORTEVENT;
0188     }
0189 
0190   if(std::isnan(zvtx))
0191     {
0192       if(_debug > 0) cout << "zvtx is NAN after attempting to grab it. ABORT EVENT!" << endl;
0193       return Fun4AllReturnCodes::ABORTEVENT;
0194     }
0195   if(abs(zvtx) > 300)
0196     {
0197       if(_debug > 2) cout << "zvtx > 300. zvtx: " << zvtx << endl;
0198       return Fun4AllReturnCodes::ABORTEVENT;
0199     }
0200 
0201   _mbzhist->Fill(zvtx);
0202 
0203   return Fun4AllReturnCodes::EVENT_OK;
0204     
0205 }
0206 //____________________________________________________________________________..
0207 int triggercountmodule::ResetEvent(PHCompositeNode *topNode)
0208 {
0209   if (Verbosity() > 0)
0210     {
0211       std::cout << "triggercountmodule::ResetEvent(PHCompositeNode *topNode) Resetting internal structures, prepare for next event" << std::endl;
0212     }
0213   return Fun4AllReturnCodes::EVENT_OK;
0214 }
0215 
0216 //____________________________________________________________________________..
0217 int triggercountmodule::EndRun(const int runnumber)
0218 {
0219   if (Verbosity() > 0)
0220     {
0221       std::cout << "triggercountmodule::EndRun(const int runnumber) Ending Run for Run " << runnumber << std::endl;
0222     }
0223   return Fun4AllReturnCodes::EVENT_OK;
0224 }
0225 
0226 //____________________________________________________________________________..
0227 int triggercountmodule::End(PHCompositeNode *topNode)
0228 {
0229   if (Verbosity() > 0)
0230     {
0231       std::cout << "triggercountmodule::End(PHCompositeNode *topNode) This is the End..." << std::endl;
0232     }
0233   //cout << _evtn << endl;
0234   if(_nseg != _lastseg && _evtn < 49500)
0235     {
0236       _badFlag = 1;
0237     }
0238 
0239   for(int i=0; i<64; ++i)
0240     {
0241       if(_endScal[i] == _startScal[i]) _avgPS[i] = 0;
0242       else _avgPS[i] = ((float)(_endLive[i]-_startLive[i]))/(_endScal[i]-_startScal[i]);
0243       if(_debug > 2) cout << "realend" << _endLive[i] << endl;
0244     }
0245   for(int i=0; i<5; ++i)
0246     {
0247       _eMBDlive[i] = _avgPS[10]*_trigCounts[i][10];
0248     }
0249 
0250   int avgBunchCounts = _bunchHist->Integral() / 121.0;
0251   
0252   for(int i=1; i<122; ++i)
0253     {
0254       if(_bunchHist->GetBinContent(i) > avgBunchCounts/5)
0255     {
0256       _nBunch++;
0257     }
0258     }
0259   _outfile = new TFile(_filename.c_str(), "RECREATE");
0260   _outfile->cd();
0261   _bunchHist->Write();
0262   _mbzhist->Write();
0263   if(_debug > 0) _tree->Print();
0264   _tree->Fill();
0265   _tree->Write();
0266   _outfile->Write();
0267   _outfile->Close();
0268 
0269 
0270   return Fun4AllReturnCodes::EVENT_OK;
0271 }
0272 
0273 //____________________________________________________________________________..
0274 int triggercountmodule::Reset(PHCompositeNode *topNode)
0275 {
0276   if (Verbosity() > 0)
0277     {
0278       std::cout << "triggercountmodule::Reset(PHCompositeNode *topNode) being Reset" << std::endl;
0279     }
0280   return Fun4AllReturnCodes::EVENT_OK;
0281 }
0282 
0283 //____________________________________________________________________________..
0284 void triggercountmodule::Print(const std::string &what) const
0285 {
0286   std::cout << "triggercountmodule::Print(const std::string &what) const Printing info for " << what << std::endl;
0287 }