File indexing completed on 2025-08-14 08:14:26
0001 #include "laserQA.h"
0002
0003 #include <trackbase/LaserClusterContainer.h>
0004 #include <trackbase/LaserCluster.h>
0005 #include <trackbase/TpcDefs.h>
0006
0007
0008 #include <ffaobjects/EventHeader.h>
0009
0010 #include <fun4all/Fun4AllReturnCodes.h>
0011
0012 #include <phool/PHCompositeNode.h>
0013 #include <phool/getClass.h>
0014 #include <phool/phool.h>
0015
0016
0017
0018 #include <TVector3.h>
0019
0020 #include <string>
0021
0022
0023
0024 laserQA::laserQA(const std::string& name)
0025 : SubsysReco(name)
0026 {
0027 }
0028
0029 int laserQA::InitRun(PHCompositeNode* topNode)
0030 {
0031
0032 outputFile = new TFile(m_output.c_str(),"RECREATE");
0033
0034
0035
0036 T = new TTree("T","T");
0037 T->Branch("hits",&hitHist);
0038 T->Branch("adc",&adc);
0039 T->Branch("layer",&layer);
0040 T->Branch("iphi",&iphi);
0041 T->Branch("it",&it);
0042 T->Branch("phi",&phi);
0043 T->Branch("R",&R);
0044
0045
0046
0047 return Fun4AllReturnCodes::EVENT_OK;
0048
0049 }
0050
0051
0052 int laserQA::process_event(PHCompositeNode* topNode)
0053 {
0054
0055 LaserClusterContainer *lcc = findNode::getClass<LaserClusterContainer>(topNode, "LASER_CLUSTER");
0056 if (!lcc)
0057 {
0058 std::cout << PHWHERE << "LASER_CLUSTER Node missing, abort." << std::endl;
0059 return Fun4AllReturnCodes::ABORTRUN;
0060 }
0061
0062
0063 if(lcc->size() < 1000) return Fun4AllReturnCodes::EVENT_OK;
0064
0065
0066
0067
0068
0069
0070
0071
0072 int clusNum = -1;
0073
0074 auto clusrange = lcc->getClusters();
0075 for (auto cmitr = clusrange.first;
0076 cmitr != clusrange.second;
0077 ++cmitr)
0078 {
0079
0080
0081
0082
0083 clusNum++;
0084
0085
0086 const auto& [cmkey, cmclus_orig] = *cmitr;
0087 LaserCluster* cmclus = cmclus_orig;
0088
0089
0090
0091 layer = cmclus->getLayer();
0092
0093
0094
0095
0096 adc = cmclus->getAdc();
0097 iphi = cmclus->getIPhi();
0098 it = cmclus->getIT();
0099
0100 float X = cmclus->getX();
0101 float Y = cmclus->getY();
0102
0103 R = sqrt(X*X + Y*Y);
0104 phi = atan2(Y,X);
0105
0106 hitHist = new TH3D("hitHist",";layer;iphi;it;adc",9,round(layer)-4.5,round(layer)+4.5,9,round(iphi)-4.5,round(iphi)+4.5,15,round(it)-7.5,round(it)+7.5);
0107
0108
0109 unsigned int nHits = cmclus->getNhits();
0110
0111 for(unsigned int i=0; i<nHits; i++)
0112 {
0113
0114 hitHist->Fill(cmclus->getHitLayer(i), cmclus->getHitIPhi(i), cmclus->getHitIT(i), cmclus->getHitAdc(i));
0115 }
0116
0117
0118 T->Fill();
0119
0120 if(hitHist)
0121 {
0122 delete hitHist;
0123 }
0124
0125 }
0126
0127 return Fun4AllReturnCodes::EVENT_OK;
0128 }
0129
0130 int laserQA::End(PHCompositeNode* )
0131 {
0132 outputFile->cd();
0133
0134 T->Write();
0135
0136 outputFile->Close();
0137
0138 return Fun4AllReturnCodes::EVENT_OK;
0139 }