Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:15:02

0001 /** @file PHAjMaker.h
0002     @author Kolja Kauder
0003     @version Revision 0.2
0004     @brief Class for A<SUB>J</SUB> analysis
0005     @details Uses JetAnalyzer objects to perform A<SUB>J</SUB> analysis.
0006     @date Aug 03, 2015
0007 */
0008 
0009 #ifndef PHAJMAKER_H__
0010 #define PHAJMAKER_H__
0011 
0012 #include <fun4all/SubsysReco.h>
0013 #include <fun4all/Fun4AllServer.h>
0014 #include <fun4all/Fun4AllHistoManager.h>
0015 #include <map>
0016 #include <set>
0017 #include <string>
0018 #include <vector>
0019 
0020 #include "PJTranslator.h"
0021 
0022 #ifndef __CINT__
0023 #include "JetAnalyzer.hh"
0024 #include "AjParameters.hh"
0025 #endif // __CINT__
0026 
0027 #include "TH1.h"
0028 #include "TH2.h"
0029 #include "TH3.h"
0030 #include "TString.h"
0031 #include "TChain.h"
0032 
0033 class PHCompositeNode;
0034 class RawClusterContainer;
0035 class SvtxTrackMap;
0036 class TF1;
0037 
0038 #include <assert.h>
0039 #include <iostream>
0040 #include <cmath>
0041 
0042 /* namespace fastjet { */
0043 /*   class PseudoJet; */
0044 /*   class JetDefinition; */
0045 /*   class SISConePlugin;     */
0046 /*   class ClusterSequence; */
0047 /* } */
0048 
0049 
0050 class PHAJMaker: public SubsysReco
0051 {
0052   // ------------------------------------------------------
0053   // ------------------- PUBLIC  --------------------------
0054   // ------------------------------------------------------
0055  public:
0056 
0057   /** Standard constructor. Set up analysis parameters.
0058       \param name: Name
0059       \param R: jet resolution parameter (radius)
0060       \param jet_ptmin: minimum jet p<SUB>T</SUB>
0061       \param jet_ptmax: maximum jet p<SUB>T</SUB>
0062       \param LeadPtMin: leading jet minimum p<SUB>T</SUB>
0063       \param SubLeadPtMin: subleading jet minimum p<SUB>T</SUB>
0064       \param max_track_rap: constituent rapidity cut
0065       \param PtConsLo: constituent minimum p<SUB>T</SUB>
0066       \param PtConsHi: constituent maximum p<SUB>T</SUB>
0067       \param dPhiCut: opening angle for dijet requirement. Accept only  |&phi;1 - &phi;2 - &pi;| < dPhiCut.
0068    */
0069   PHAJMaker( const std::string &name = "PHAJMaker",
0070          double R = 0.4,
0071          double jet_ptmin = 10.0, double jet_ptmax = 100.0,
0072          double LeadPtMin = 20.0, double SubLeadPtMin = 10, 
0073          double max_track_rap = 1.0, double PtConsLo=0.2, double PtConsHi=2.0,
0074          double dPhiCut = 0.4
0075          );
0076 
0077   /** Initializer, called ONCE */
0078   int Init(PHCompositeNode *);
0079 
0080   /** Main analysis routine, called for EVERY EVENT */
0081   int process_event(PHCompositeNode *);
0082 
0083   /** This little helper is true if there's at least one 10 GeV jet
0084    **/
0085   bool Has10Gev;
0086 
0087   // Getters and Setters
0088   // -------------------
0089   /// Get jet radius
0090   inline double GetR ( )                   { return R; };
0091   /// Set jet radius
0092   inline void   SetR ( const double newv ) { R=newv;   };
0093 
0094   /// Get minimum jet p<SUB>T</SUB>
0095   inline double GetJet_ptmin ( )                   { return jet_ptmin; };
0096   /// Set minimum jet p<SUB>T</SUB>
0097   inline void   SetJet_ptmin ( const double newv ) { jet_ptmin=newv;   };
0098 
0099   /// Get maximum jet p<SUB>T</SUB>
0100   inline double GetJet_ptmax ( )                   { return jet_ptmax; };
0101   /// Set maximum jet p<SUB>T</SUB>
0102   inline void   SetJet_ptmax ( const double newv ) { jet_ptmax=newv;   };
0103 
0104   /// Get jet rapidity acceptance
0105   inline double GetMax_rap ( )                   { return max_rap; };
0106   /// Set jet rapidity acceptance
0107   inline void   SetMax_rap ( const double newv ) { max_rap=newv;   };
0108 
0109   /// Get ghosted area rapidity cut, should be >= max_rap + 2*R
0110   inline double GetGhost_maxrap ( )                   { return ghost_maxrap; };
0111   /// Set ghosted area rapidity cut, should be >= max_rap + 2*R
0112   inline void   SetGhost_maxrap ( const double newv ) { ghost_maxrap=newv;   };
0113 
0114    /// Get dijet opening angle
0115   inline double GetDPhiCut ( )                   { return dPhiCut; };
0116   /// Set dijet opening angle
0117   inline void   SetDPhiCut ( const double newv ) { dPhiCut=newv;   };  
0118 
0119   /** Standard dtor (empty) */
0120   virtual ~PHAJMaker();
0121 
0122   /** Dijet asymmetry A<SUB>J</SUB> = &Delta;p<SUB>T</SUB> / &Sigma;p<SUB>T</SUB> */
0123 
0124 #ifndef __CINT__
0125   inline double CalcAj ( std::vector<fastjet::PseudoJet>& jets ){
0126     if ( jets.size()!=2 ){
0127       throw ( -1 );
0128       return -1e10;
0129     }
0130     return fabs (( jets.at(0).pt()-jets.at(1).pt() ) / ( jets.at(0).pt()+jets.at(1).pt() ));    
0131   }
0132 
0133   // Objects will be handed by _reference_! Obviates need for setter
0134   /// Handle to jet definition
0135   inline fastjet::JetDefinition& GetJet_def () { return jet_def; }
0136   /// Handle to selector for low  p<SUB>T</SUB> constituents
0137   inline fastjet::Selector& GetLoConsSelector () { return slo; }
0138   /// Handle to selector for high  p<SUB>T</SUB> constituents
0139   inline fastjet::Selector& GetHiConsSelector () { return shi; }
0140   
0141   /// Handle to selector for jet candidates
0142   inline fastjet::Selector& GetJetSelector () { return sjet; }
0143 
0144   /// Handle to ghosted area specification
0145   inline fastjet::GhostedAreaSpec& GetArea_spec () { return area_spec; }
0146   /// Handle to jet area definition
0147   inline fastjet::AreaDefinition& GetArea_def () { return area_def; }
0148 
0149 
0150   /// Handle to JetAnalyzer for high pT constituents
0151   inline JetAnalyzer* GetJAhi() {return pJAhi; };
0152   /// Handle to JetAnalyzer for low pT constituents
0153   inline JetAnalyzer* GetJAlo() {return pJAlo; };
0154 
0155   /// Handle to unaltered clustering result with high pT constituents
0156   inline std::vector<fastjet::PseudoJet> GetJAhiResult() {return JAhiResult; };
0157   /// Handle to unaltered clustering result with low pT constituents
0158   inline std::vector<fastjet::PseudoJet> GetJAloResult() {return JAloResult; };
0159 
0160   /// Handle to high pT constituents
0161   inline std::vector<fastjet::PseudoJet> GetLoConstituents() {return pLo; };
0162   /// Handle to low pT constituents
0163   inline std::vector<fastjet::PseudoJet> GetHiConstituents() {return pHi; };
0164 
0165   
0166   /// Handle to Dijet result with high pT constituents
0167   inline std::vector<fastjet::PseudoJet> GetDiJetsHi() {return DiJetsHi; };
0168   /// Handle to Dijet result with low pT constituents
0169   inline std::vector<fastjet::PseudoJet> GetDiJetsLo() {return DiJetsLo; };
0170 
0171 
0172 #endif // __CINT__
0173 
0174   //void save_jets_to_nodetree();
0175 
0176 
0177   
0178   // ------------------------------------------------------
0179   // ------------------- PRIVATE --------------------------
0180   // ------------------------------------------------------
0181   private:           
0182 
0183   Fun4AllHistoManager* MyHistos;
0184   TH2D* UnmatchedAJ_hi;
0185   TH2D* AJ_hi;
0186   TH2D* AJ_lo;
0187 
0188   double R;              ///< Resolution parameter ("jet radius")
0189   double jet_ptmin;      ///< minimum jet p<SUB>T</SUB>
0190   double jet_ptmax;      ///< maximum jet p<SUB>T</SUB>
0191   double LeadPtMin;      ///< leading jet minimum p<SUB>T</SUB>
0192   double SubLeadPtMin;   ///< subleading jet minimum p<SUB>T</SUB>
0193   double max_rap;        ///< jet rapidity acceptance
0194   double ghost_maxrap;   ///< for ghosted area, should be >= max_rap + 2*R
0195 
0196   double max_track_rap;  ///< constituent rapidity cut
0197   double PtConsLo;       ///< constituent minimum p<SUB>T</SUB>
0198   double PtConsHi;       ///< constituent maximum p<SUB>T</SUB>
0199 
0200   double dPhiCut;        ///< opening angle for dijet requirement. Accept only  |&phi;1 - &phi;2 - &pi;| < &Delta;&phi;.
0201 
0202 #ifndef __CINT__
0203   fastjet::JetDefinition jet_def;       ///< jet definition
0204   fastjet::JetDefinition other_jet_def; ///< jet definition with a different radius
0205 
0206   fastjet::Selector select_track_rap;   ///< constituent rapidity selector
0207   fastjet::Selector select_lopt;        ///< constituent p<SUB>T</SUB> selector
0208   fastjet::Selector select_hipt;        ///< constituent p<SUB>T</SUB> selector
0209 
0210   fastjet::Selector slo;                ///< compound selector for low  p<SUB>T</SUB> constituents
0211   fastjet::Selector shi;                ///< compound selector for high p<SUB>T</SUB> constituents
0212 
0213 // Relevant jet candidates
0214   fastjet::Selector select_jet_rap;        ///< jet rapidity selector
0215   fastjet::Selector select_jet_pt_min;     ///< jet p<SUB>T</SUB> selector
0216   fastjet::Selector select_jet_pt_max;     ///< jet p<SUB>T</SUB> selector
0217   fastjet::Selector sjet;                  ///< compound jet selector
0218 
0219   fastjet::GhostedAreaSpec area_spec;      ///< ghosted area specification
0220   fastjet::AreaDefinition area_def;        ///< jet area definition
0221 
0222   JetAnalyzer* pJAhi;                      ///< JetAnalyzer object for high pT
0223   JetAnalyzer* pJAlo;                      ///< JetAnalyzer object for low pT
0224   JetAnalyzer* pOtherJAlo;                 ///< JetAnalyzer object for low pT with different R
0225 
0226   std::vector<fastjet::PseudoJet> pHi;     ///< High pT constituents
0227   std::vector<fastjet::PseudoJet> pLo;     ///< Low pT constituents
0228 
0229   std::vector<fastjet::PseudoJet> JAhiResult;  ///< Unaltered clustering result with high pT constituents
0230   std::vector<fastjet::PseudoJet> JAloResult;  ///< Unaltered clustering result with low pT constituents
0231   std::vector<fastjet::PseudoJet> OtherJAloResult;  ///< Unaltered clustering result with low pT constituents, different R
0232   
0233   std::vector<fastjet::PseudoJet> DiJetsHi;    ///< Dijet result with high pT constituents
0234   std::vector<fastjet::PseudoJet> DiJetsLo;    ///< Dijet result with low pT constituents
0235   std::vector<fastjet::PseudoJet> OtherDiJetsLo; ///< Dijet result with low pT constituents and different R
0236 
0237 #endif // __CINT__
0238 
0239 
0240 };
0241 
0242 #endif