Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:16:11

0001 /*
0002  * TPCFEETestRecov1.h
0003  *
0004  *  Created on: Sep 19, 2018
0005  *      Author: jinhuang
0006  */
0007 
0008 #ifndef CORESOFTWARE_OFFLINE_PACKAGES_TPCDAQ_TPCFEETESTRECOV1_H_
0009 #define CORESOFTWARE_OFFLINE_PACKAGES_TPCDAQ_TPCFEETESTRECOV1_H_
0010 
0011 #include <fun4all/SubsysReco.h>
0012 
0013 #include <TObject.h>
0014 
0015 #include <stdint.h>
0016 #include <cmath>
0017 #include <map>
0018 #include <set>
0019 #include <string>
0020 #include <vector>
0021 
0022 class PHCompositeNode;
0023 class Fun4AllHistoManager;
0024 class TTree;
0025 class TClonesArray;
0026 class Event;
0027 
0028 namespace TPCDaqDefs
0029 {
0030 namespace FEEv1
0031 {
0032 class SampleFit_PowerLawDoubleExp_PDFMaker;
0033 }
0034 }  // namespace TPCDaqDefs
0035 
0036 class TPCFEETestRecov1 : public SubsysReco
0037 {
0038  public:
0039   TPCFEETestRecov1(const std::string &outputfilename =
0040                        "TPCFEETestRecov1.root");
0041   virtual ~TPCFEETestRecov1();
0042 
0043   int Init(PHCompositeNode *topNode);
0044   int InitRun(PHCompositeNode *topNode);
0045   int process_event(PHCompositeNode *topNode);
0046   int ResetEvent(PHCompositeNode *topNode);
0047   int End(PHCompositeNode *topNode);
0048 
0049   void setClusteringZeroSuppression(int threshold)
0050   {
0051     m_clusteringZeroSuppression = threshold;
0052   }
0053 
0054   void setNPostSample(int nPostSample)
0055   {
0056     m_nPostSample = nPostSample;
0057   }
0058 
0059   void setNPreSample(int nPreSample)
0060   {
0061     m_nPreSample = nPreSample;
0062   }
0063 
0064   //! simple event header class for ROOT file IO
0065   class EventHeader : public TObject
0066   {
0067    public:
0068     int run;
0069     int event;
0070 
0071     uint32_t bx_counter;
0072     bool bx_counter_consistent;
0073 
0074     int xray_x;
0075     int xray_y;
0076 
0077     EventHeader()
0078       : run(-1)
0079       , event(-1)
0080       , bx_counter(0)
0081       , bx_counter_consistent(true)
0082       , xray_x(-1)
0083       , xray_y(-1)
0084     {
0085     }
0086 
0087     ClassDefOverride(TPCFEETestRecov1::EventHeader, 1)
0088   };
0089 
0090   //! buffer for full event data
0091   class PadPlaneData
0092   {
0093    public:
0094     PadPlaneData();
0095     void Reset();
0096 
0097     struct SampleID
0098     {
0099       int pady;
0100       int padx;
0101       int sample;
0102 
0103       void adjust(const SampleID &adjustment)
0104       {
0105         pady += adjustment.pady;
0106         padx += adjustment.padx;
0107         sample += adjustment.sample;
0108       }
0109     };
0110 
0111     static bool IsValidPad(const int pad_x, const int pad_y);
0112     std::vector<int> &getPad(const int pad_x, const int pad_y);
0113     int getSample(const SampleID &id);
0114 
0115     //! 3-D Graph clustering based on PHMakeGroups()
0116     void Clustering(int zero_suppression, bool verbosity = false);
0117 
0118 #if !defined(__CINT__) || defined(__CLING__)
0119 
0120     const std::vector<std::vector<std::vector<int>>> &getData() const
0121     {
0122       return m_data;
0123     }
0124 
0125     const std::multimap<int, SampleID> &getGroups() const
0126     {
0127       return m_groups;
0128     }
0129 
0130    private:
0131     //! full event data in index order of m_data[pady][padx][sample]
0132     std::vector<std::vector<std::vector<int>>> m_data;
0133 
0134     std::multimap<int, SampleID> m_groups;
0135 
0136 #endif  // #if !defined(__CINT__) || defined(__CLING__)
0137   };
0138 
0139   //! buffer for a cluster's data
0140   class ClusterData : public TObject
0141   {
0142    public:
0143     ClusterData()
0144       : min_sample(-1)
0145       , max_sample(-1)
0146       , peak(NAN)
0147       , peak_sample(NAN)
0148       , pedstal(NAN)
0149       , avg_padx(NAN)
0150       , avg_pady(NAN)
0151       , size_pad_x(-1)
0152       , size_pad_y(-1)
0153     {
0154     }
0155 
0156     std::set<int> padxs;
0157     std::set<int> padys;
0158     std::set<int> samples;
0159 
0160 #if !defined(__CINT__) || defined(__CLING__)
0161     std::map<int, std::vector<double>> padx_samples;
0162     std::map<int, std::vector<double>> pady_samples;
0163     std::vector<double> sum_samples;
0164 #endif  // #if !defined(__CINT__) || defined(__CLING__)
0165 
0166     int min_sample;
0167     int max_sample;
0168 
0169     double peak;
0170     double peak_sample;
0171     double pedstal;
0172 
0173     std::map<int, double> padx_peaks;
0174     std::map<int, double> pady_peaks;
0175 
0176     //! pad coordinate
0177     double avg_padx;
0178     double avg_pady;
0179 
0180     //! pad size
0181     int size_pad_x;
0182     int size_pad_y;
0183 
0184     ClassDefOverride(ClusterData, 1);
0185   };
0186 
0187   //! simple channel header class for ROOT file IO
0188   class ChannelHeader : public TObject
0189   {
0190    public:
0191     int size;
0192     //! = p->iValue(channel * kPACKET_LENGTH + 2) & 0xffff;  // that's the Elink packet type
0193     uint8_t packet_type;
0194     //! = ((p->iValue(channel * kPACKET_LENGTH + 4) & 0xffff) << 4) | (p->iValue(channel * kPACKET_LENGTH + 5) & 0xffff);
0195     uint32_t bx_counter;
0196     //! = (p->iValue(channel * kPACKET_LENGTH + 3) >> 5) & 0xf;
0197     uint8_t sampa_address;
0198     //! = p->iValue(channel * kPACKET_LENGTH + 3) & 0x1f;
0199     uint16_t sampa_channel;
0200     //! = (sampa_address << 5) | sampa_channel;
0201     uint16_t fee_channel;
0202 
0203     //! pad coordinate
0204     int pad_x;
0205     int pad_y;
0206 
0207     int pedestal;
0208     int max;
0209 
0210     ChannelHeader()
0211       : size(0)
0212       , packet_type(0)
0213       , bx_counter(0)
0214       , sampa_address(0)
0215       , sampa_channel(0)
0216       , fee_channel(0)
0217       , pad_x(-1)
0218       , pad_y(-1)
0219       , pedestal(-1)
0220       , max(-1)
0221     {
0222     }
0223 
0224     ClassDefOverride(TPCFEETestRecov1::ChannelHeader, 1)
0225   };
0226 
0227  private:
0228 #if !defined(__CINT__) || defined(__CLING__)
0229 
0230   // IO stuff
0231 
0232   Fun4AllHistoManager *getHistoManager();
0233 
0234   std::string m_outputFileName;
0235 
0236   TTree *m_eventT;
0237 
0238   EventHeader m_eventHeader;
0239   EventHeader *m_peventHeader;  //! ->m_eventHeader,  for filling TTree
0240 
0241   int m_nClusters;
0242   TClonesArray *m_IOClusters;
0243 
0244   TTree *m_chanT;
0245 
0246   ChannelHeader m_chanHeader;
0247   ChannelHeader *m_pchanHeader;  //! ->m_chanHeader,  for filling TTree
0248   std::vector<uint32_t> m_chanData;
0249 
0250   // clustering stuff
0251   PadPlaneData m_padPlaneData;
0252   std::map<int, ClusterData> m_clusters;
0253 
0254   //! rough zero suppression by subtracting sample medium value
0255   //! \return pair of pedestal and max
0256   static std::pair<int, int> roughZeroSuppression(std::vector<int> &data);
0257 
0258   //! Clustering then prepare IOs
0259   void Clustering(void);
0260 
0261 #endif  // #if !defined(__CINT__) || defined(__CLING__)
0262 
0263   int m_clusteringZeroSuppression;
0264   int m_nPreSample;
0265   int m_nPostSample;
0266 
0267   void get_motor_loc(Event *evt);
0268 
0269   int m_XRayLocationX;
0270   int m_XRayLocationY;
0271 
0272   TPCDaqDefs::FEEv1::SampleFit_PowerLawDoubleExp_PDFMaker *m_pdfMaker;
0273 };
0274 
0275 bool operator<(const TPCFEETestRecov1::PadPlaneData::SampleID &s1, const TPCFEETestRecov1::PadPlaneData::SampleID &s2);
0276 
0277 #endif /* CORESOFTWARE_OFFLINE_PACKAGES_TPCDAQ_TPCFEETESTRECOV1_H_ */