|
||||
File indexing completed on 2024-09-07 04:37:17
0001 #ifndef PhysicsTools_Heppy_Hemisphere_h 0002 #define PhysicsTools_Heppy_Hemisphere_h 0003 0004 /* \class Hemisphere 0005 * 0006 * Class that, given the 4-momenta of the objects in the event, 0007 * allows to split up the event into two "hemispheres" according 0008 * to different criteria (see below)/ 0009 * 0010 * Authors: Luc Pape & Filip Moortgat Date: July 2005 0011 * Updated: July 2006 0012 * Updated: 15 Sept 2006 0013 * Updated: 15 November 2006 0014 * Updated: 01 November 2008 0015 * Updated: 13 December 2010 (P.Nef) 0016 * Updated: 09 February 2011 0017 * Updated: 26 September 2011 (P.Nef) 0018 */ 0019 0020 #include <vector> 0021 #include <iostream> 0022 #include <cmath> 0023 0024 namespace heppy { 0025 0026 class Hemisphere { 0027 public: 0028 // There are 2 constructors: 0029 // 1. Constructor taking as argument vectors of Px, Py, Pz and E of the objects in 0030 // the event that should be separated, the seeding method and the hemisphere 0031 // association method, 0032 // 2. Constructor taking as argument vectors of Px, Py, Pz and E of the objects in 0033 // the event that should be separated. The seeding method and the hemisphere 0034 // association method should then be defined by SetMethod(seeding_method, association_method). 0035 // 0036 // Seeding method: choice of 2 inital axes 0037 // 1: 1st: max P ; 2nd: max P * delta R wrt first one (if delta R > 0.5) 0038 // 2: 2 objects who give maximal invariant mass (recommended) 0039 // 3: 2 objects who give maximal transverse mass 0040 // 4: 2 leading jets (min dR can be set by SetDRminSeed1()) 0041 // 0042 // Hemisphere association method: 0043 // 1: maximum pt longitudinal projected on the axes 0044 // 2: minimal mass squared sum of the hemispheres 0045 // 3: minimal Lund distance (recommended) 0046 // 0047 // Note that SetMethod() also allows the seeding and/or association method to be 0048 // redefined for an existing hemisphere object. The GetAxis() or GetGrouping() is 0049 // then recomputed using the newly defined methods. 0050 // 0051 // Some parameters possibly affecting the logic of hemisphere reconstruction can also be set 0052 // by the user: 0053 // SetNoSeed() to prevent the given object from being used as a seed 0054 // (but it can be associated to hemispheres). 0055 // SetNoAssoc() to prevent the given object from being used in the hemisphere association 0056 // (then it cannot be used as seed either). 0057 // ClearAllNoLists() to reset the list of NoSeed and NoAssoc objects to empty. 0058 // SetDRminSeed1() minimum value of DeltaR between the two hemisphere seed for seed method 1 0059 // (default = 0.5). Does not affect other methods. 0060 // SetnItermax() maximum number of iterations allowed for the association methods (default = 100). 0061 // 0062 // Methods are provided to avoid including ISR jets into the hemispheres. 0063 // This is done by declaring the variable type and cut value above which an object should not 0064 // be included in the hemispheres. 0065 // Only 1 cut can be used at the time (a new declaration overwrites the previous one). 0066 // RejectISRPtmax() max Pt w.r.t. the hemisphere below which objects can be included 0067 // (default = 10000. GeV) 0068 // RejectISRDRmax() max DeltaR below which objects can be included 0069 // (default = 10.) 0070 0071 Hemisphere() {} 0072 0073 Hemisphere(std::vector<float> Px_vector, 0074 std::vector<float> Py_vector, 0075 std::vector<float> Pz_vector, 0076 std::vector<float> E_vector, 0077 int seed_method, 0078 int hemisphere_association_method); 0079 0080 Hemisphere(std::vector<float> Px_vector, 0081 std::vector<float> Py_vector, 0082 std::vector<float> Pz_vector, 0083 std::vector<float> E_vector); 0084 0085 // Destructor 0086 ~Hemisphere() {} 0087 0088 // return Nx, Ny, Nz, P, E of the axis of group 1 0089 std::vector<float> getAxis1(); 0090 // return Nx, Ny, Nz, P, E of the axis of group 2 0091 std::vector<float> getAxis2(); 0092 0093 // where Nx, Ny, Nz are the direction cosines e.g. Nx = Px/P, 0094 // P is the momentum, E is the energy 0095 0096 // return vector with "1" and "2"'s according to which group the object belongs 0097 // and 0 if the object is not associated to any hemisphere 0098 // (order of objects in vector is same as input) 0099 std::vector<int> getGrouping(); 0100 0101 // set or overwrite the seed and association methods 0102 void SetMethod(int seed_method, int hemisphere_association_method) { 0103 seed_meth = seed_method; 0104 hemi_meth = hemisphere_association_method; 0105 status = 0; 0106 } 0107 0108 // prevent an object from being used as a seed (but it can be associated to hemispheres) 0109 // (method introduced on 15/09/06) 0110 void SetNoSeed(int object_number) { 0111 Object_Noseed[object_number] = 1; 0112 status = 0; 0113 } 0114 0115 // prevent an object from being used for hemisphere asociation (and for seeding) 0116 // (method introduced on 01/11/08) 0117 void SetNoAssoc(int object_number) { 0118 Object_Noassoc[object_number] = 1; 0119 Object_Noseed[object_number] = 1; 0120 status = 0; 0121 } 0122 0123 // reset the list of NoSeed and NoAssoc objects to empty 0124 // (method introduced on 01/11/08) 0125 void ClearAllNoLists() { 0126 for (int i = 0; i < (int)Object_Noseed.size(); ++i) { 0127 Object_Noassoc[i] = 0; 0128 Object_Noseed[i] = 0; 0129 status = 0; 0130 } 0131 } 0132 0133 // set the min DeltaR for the seeds of seed method 1 0134 // (method introduced on 01/11/08) 0135 void SetDRminSeed1(float rmin) { dRminSeed1 = rmin; } 0136 0137 // set the maximum number of iterations allowed for association 0138 // (method introduced on 01/11/08) 0139 void SetnItermax(int niter) { nItermax = niter; } 0140 0141 // set max Pt w.r.t. the hemisphere below which objects can be included 0142 void RejectISRPtmax(float ptmax) { 0143 rejectISRPtmax = ptmax; 0144 rejectISRPt = 1; 0145 rejectISRDR = 0; 0146 rejectISR = 1; 0147 status = 0; 0148 } 0149 0150 // set max DeltaR below which objects can be included 0151 void RejectISRDRmax(float drmax) { 0152 rejectISRDRmax = drmax; 0153 rejectISRDR = 1; 0154 rejectISRPt = 0; 0155 rejectISR = 1; 0156 status = 0; 0157 } 0158 0159 // controls the level of debug prints 0160 void SetDebug(int debug) { dbg = debug; } 0161 int GetNumLoop() { return numLoop; } 0162 0163 private: 0164 // the hemisphere separation algorithm 0165 int Reconstruct(); 0166 int RejectISR(); 0167 0168 std::vector<float> Object_Px; 0169 std::vector<float> Object_Py; 0170 std::vector<float> Object_Pz; 0171 std::vector<float> Object_P; 0172 std::vector<float> Object_Pt; 0173 std::vector<float> Object_E; 0174 std::vector<float> Object_Phi; 0175 std::vector<float> Object_Eta; 0176 std::vector<int> Object_Group; 0177 std::vector<int> Object_Noseed; 0178 std::vector<int> Object_Noassoc; 0179 0180 std::vector<float> Axis1; 0181 std::vector<float> Axis2; 0182 0183 //static const float hemivsn = 1.01; 0184 int seed_meth; 0185 int hemi_meth; 0186 int status; 0187 float dRminSeed1; 0188 int nItermax; 0189 int rejectISR; 0190 int rejectISRPt; 0191 float rejectISRPtmax; 0192 int rejectISRDR; 0193 float rejectISRDRmax; 0194 int dbg; 0195 int numLoop; 0196 }; 0197 } // namespace heppy 0198 0199 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.2.1 LXR engine. The LXR team |