Go to the documentation of this file.00001 #include "PhysicsTools/PatAlgos/interface/HemisphereAlgo.h"
00002
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004
00005 #include "DataFormats/Math/interface/deltaPhi.h"
00006 #include "DataFormats/Math/interface/deltaR.h"
00007
00008 using namespace std;
00009 using std::vector;
00010
00011 HemisphereAlgo::HemisphereAlgo(const std::vector<reco::CandidatePtr>& componentPtr_,const int seed_method, const int hemisphere_association_method )
00012 : Object(componentPtr_), Object_Group() , Axis1(), Axis2(), seed_meth(seed_method), hemi_meth(hemisphere_association_method), status(0) {
00013
00014 for(int i = 0; i < (int) Object.size(); i++){
00015 Object_Noseed.push_back(0);
00016 }
00017
00018 }
00019
00020
00021
00022 vector<float> HemisphereAlgo::getAxis1(){
00023 if (status != 1){this->reconstruct();}
00024 return Axis1;
00025 }
00026 vector<float> HemisphereAlgo::getAxis2(){
00027 if (status != 1){this->reconstruct();}
00028 return Axis2;
00029 }
00030
00031 vector<int> HemisphereAlgo::getGrouping(){
00032 if (status != 1){this->reconstruct();}
00033 return Object_Group;
00034 }
00035
00036 int HemisphereAlgo::reconstruct(){
00037
00038 int vsize = (int) Object.size();
00039
00040 LogDebug("HemisphereAlgo") << " HemisphereAlgo method ";
00041
00042 Object_Group.clear();
00043 Axis1.clear();
00044 Axis2.clear();
00045
00046 for(int j = 0; j < vsize; j++){
00047 Object_Group.push_back(0);
00048 }
00049
00050 for(int j = 0; j < 5; j++){
00051 Axis1.push_back(0);
00052 Axis2.push_back(0);
00053 }
00054
00055
00056 for (int i = 0; i <vsize; i++){
00057
00058 if ( (*(Object)[i]).p() > (*Object[i]).energy() + 0.001) {
00059
00060 edm::LogWarning("HemisphereAlgo") << "Object " << i << " has E = " << (*Object[i]).energy()
00061 << " less than P = " << (*Object[i]).p() ;
00062
00063 }
00064 }
00065
00066
00067
00068 LogDebug("HemisphereAlgo") << " Seeding method = " << seed_meth;
00069 int I_Max = -1;
00070 int J_Max = -1;
00071
00072 if (seed_meth == 1) {
00073
00074 float P_Max = 0.;
00075 float DeltaRP_Max = 0.;
00076
00077
00078 for (int i=0;i<vsize;i++){
00079 Object_Group[i] = 0;
00080 if (Object_Noseed[i] == 0 && P_Max < (*Object[i]).p()){
00081 P_Max = (*Object[i]).p();
00082 I_Max = i;
00083 }
00084 }
00085
00086 Axis1[0] = (*Object[I_Max]).px() / (*Object[I_Max]).p();
00087 Axis1[1] = (*Object[I_Max]).py() / (*Object[I_Max]).p();
00088 Axis1[2] = (*Object[I_Max]).pz() / (*Object[I_Max]).p();
00089 Axis1[3] = (*Object[I_Max]).p();
00090 Axis1[4] = (*Object[I_Max]).energy();
00091
00092
00093 for (int i=0;i<vsize;i++){
00094
00095 float DeltaR = deltaR((*Object[i]).eta(),(*Object[i]).phi(),(*Object[I_Max]).eta(),(*Object[I_Max]).phi()) ;
00096
00097 if (Object_Noseed[i] == 0 && DeltaR > 0.5) {
00098
00099 float DeltaRP = DeltaR * (*Object[i]).p();
00100 if (DeltaRP > DeltaRP_Max){
00101 DeltaRP_Max = DeltaRP;
00102 J_Max = i;
00103 }
00104 }
00105 }
00106
00107 if (J_Max >=0){
00108 Axis2[0] = (*Object[J_Max]).px() / (*Object[J_Max]).p();
00109 Axis2[1] =(*Object[J_Max]).py() / (*Object[J_Max]).p();
00110 Axis2[2] = (*Object[J_Max]).pz() / (*Object[J_Max]).p();
00111 Axis2[3] = (*Object[J_Max]).p();
00112 Axis2[4] = (*Object[J_Max]).energy();
00113
00114 } else {
00115 return 0;
00116 }
00117 LogDebug("HemisphereAlgo") << " Axis 1 is Object = " << I_Max;
00118 LogDebug("HemisphereAlgo") << " Axis 2 is Object = " << J_Max;
00119
00120
00121 } else if ( (seed_meth == 2) | (seed_meth == 3) ) {
00122
00123 float Mass_Max = 0.;
00124 float InvariantMass = 0.;
00125
00126
00127 for (int i=0;i<vsize;i++){
00128 Object_Group[i] = 0;
00129 if (Object_Noseed[i] == 0){
00130 for (int j=i+1;j<vsize;j++){
00131 if (Object_Noseed[j] == 0){
00132
00133
00134 if (seed_meth == 2){
00135 InvariantMass = ((*Object[i]).energy() + (*Object[j]).energy())* ((*Object[i]).energy() + (*Object[j]).energy())
00136 - ((*Object[i]).px() + (*Object[j]).px())* ((*Object[i]).px() + (*Object[j]).px())
00137 - ((*Object[i]).py() + (*Object[j]).py())* ((*Object[i]).py() + (*Object[j]).py())
00138 - ((*Object[i]).pz() + (*Object[j]).pz())* ((*Object[i]).pz() + (*Object[j]).pz()) ;
00139 }
00140
00141 else if (seed_meth == 3){
00142 float pti = sqrt((*Object[i]).px()*(*Object[i]).px() + (*Object[i]).py()*(*Object[i]).py());
00143 float ptj = sqrt((*Object[j]).px()*(*Object[j]).px() + (*Object[j]).py()*(*Object[j]).py());
00144 InvariantMass = 2. * (pti*ptj - (*Object[i]).px()*(*Object[j]).px()
00145 - (*Object[i]).py()*(*Object[j]).py() );
00146 }
00147 if ( Mass_Max < InvariantMass){
00148 Mass_Max = InvariantMass;
00149 I_Max = i;
00150 J_Max = j;
00151 }
00152 }
00153 }
00154 }
00155 }
00156
00157 if (J_Max>0) {
00158
00159 Axis1[0] = (*Object[I_Max]).px() / (*Object[I_Max]).p();
00160 Axis1[1] = (*Object[I_Max]).py() / (*Object[I_Max]).p();
00161 Axis1[2] = (*Object[I_Max]).pz() / (*Object[I_Max]).p();
00162
00163 Axis1[3] = (*Object[I_Max]).p();
00164 Axis1[4] = (*Object[I_Max]).energy();
00165
00166 Axis2[0] = (*Object[J_Max]).px() / (*Object[J_Max]).p();
00167 Axis2[1] =(*Object[J_Max]).py() / (*Object[J_Max]).p();
00168 Axis2[2] = (*Object[J_Max]).pz() / (*Object[J_Max]).p();
00169
00170 Axis2[3] = (*Object[J_Max]).p();
00171 Axis2[4] = (*Object[J_Max]).energy();
00172
00173 } else {
00174 return 0;
00175 }
00176 LogDebug("HemisphereAlgo") << " Axis 1 is Object = " << I_Max;
00177 LogDebug("HemisphereAlgo") << " Axis 2 is Object = " << J_Max;
00178
00179
00180 } else {
00181 throw cms::Exception("Configuration") << "Please give a valid seeding method!";
00182 }
00183
00184
00185
00186
00187 LogDebug("HemisphereAlgo") << " Association method = " << hemi_meth;
00188
00189 int numLoop = 0;
00190 bool I_Move = true;
00191
00192
00193 while (I_Move && (numLoop < 100)){
00194
00195 I_Move = false;
00196 numLoop++;
00197 LogDebug("HemisphereAlgo") << " Iteration = " << numLoop;
00198
00199 float Sum1_Px = 0.;
00200 float Sum1_Py = 0.;
00201 float Sum1_Pz = 0.;
00202 float Sum1_P = 0.;
00203 float Sum1_E = 0.;
00204 float Sum2_Px = 0.;
00205 float Sum2_Py = 0.;
00206 float Sum2_Pz = 0.;
00207 float Sum2_P = 0.;
00208 float Sum2_E = 0.;
00209
00210
00211 if (hemi_meth == 1) {
00212
00213 for (int i=0;i<vsize;i++){
00214 float P_Long1 = (*Object[i]).px()*Axis1[0] + (*Object[i]).py()*Axis1[1] + (*Object[i]).pz()*Axis1[2];
00215 float P_Long2 = (*Object[i]).px()*Axis2[0]+ (*Object[i]).py()*Axis2[1] + (*Object[i]).pz()*Axis2[2];
00216 if (P_Long1 >= P_Long2){
00217 if (Object_Group[i] != 1){
00218 I_Move = true;
00219 }
00220 Object_Group[i] = 1;
00221 Sum1_Px += (*Object[i]).px();
00222 Sum1_Py += (*Object[i]).py();
00223 Sum1_Pz += (*Object[i]).pz();
00224 Sum1_P += (*Object[i]).p();
00225 Sum1_E += (*Object[i]).energy();
00226 } else {
00227 if (Object_Group[i] != 2){
00228 I_Move = true;
00229 }
00230 Object_Group[i] = 2;
00231 Sum2_Px += (*Object[i]).px();
00232 Sum2_Py += (*Object[i]).py();
00233 Sum2_Pz += (*Object[i]).pz();
00234 Sum2_P += (*Object[i]).p();
00235 Sum2_E += (*Object[i]).energy();
00236 }
00237 }
00238
00239 } else if (hemi_meth == 2 || hemi_meth == 3){
00240
00241 for (int i=0;i<vsize;i++){
00242 if (i == I_Max) {
00243 Object_Group[i] = 1;
00244 Sum1_Px += (*Object[i]).px();
00245 Sum1_Py += (*Object[i]).py();
00246 Sum1_Pz += (*Object[i]).pz();
00247 Sum1_P += (*Object[i]).p();
00248 Sum1_E += (*Object[i]).energy();
00249 } else if (i == J_Max) {
00250 Object_Group[i] = 2;
00251 Sum2_Px += (*Object[i]).px();
00252 Sum2_Py += (*Object[i]).py();
00253 Sum2_Pz += (*Object[i]).pz();
00254 Sum2_P += (*Object[i]).p();
00255 Sum2_E += (*Object[i]).energy();
00256 } else {
00257
00258
00259 if(!I_Move){
00260
00261 float NewAxis1_Px = Axis1[0] * Axis1[3];
00262 float NewAxis1_Py = Axis1[1] * Axis1[3];
00263 float NewAxis1_Pz = Axis1[2] * Axis1[3];
00264 float NewAxis1_E = Axis1[4];
00265
00266 float NewAxis2_Px = Axis2[0] * Axis2[3];
00267 float NewAxis2_Py = Axis2[1] * Axis2[3];
00268 float NewAxis2_Pz = Axis2[2] * Axis2[3];
00269 float NewAxis2_E = Axis2[4];
00270
00271 if (Object_Group[i] == 1){
00272
00273 NewAxis1_Px = NewAxis1_Px - (*Object[i]).px();
00274 NewAxis1_Py = NewAxis1_Py - (*Object[i]).py();
00275 NewAxis1_Pz = NewAxis1_Pz - (*Object[i]).pz();
00276 NewAxis1_E = NewAxis1_E - (*Object[i]).energy();
00277
00278 } else if (Object_Group[i] == 2) {
00279
00280 NewAxis2_Px = NewAxis2_Px - (*Object[i]).px();
00281 NewAxis2_Py = NewAxis2_Py - (*Object[i]).py();
00282 NewAxis2_Pz = NewAxis2_Pz - (*Object[i]).pz();
00283 NewAxis2_E = NewAxis2_E - (*Object[i]).energy();
00284 }
00285
00286
00287 float mass1 = NewAxis1_E - (((*Object[i]).px()*NewAxis1_Px + (*Object[i]).py()*NewAxis1_Py +
00288 (*Object[i]).pz()*NewAxis1_Pz)/(*Object[i]).p());
00289
00290 float mass2 = NewAxis2_E - (((*Object[i]).px()*NewAxis2_Px + (*Object[i]).py()*NewAxis2_Py +
00291 (*Object[i]).pz()*NewAxis2_Pz)/(*Object[i]).p());
00292
00293 if (hemi_meth == 3) {
00294
00295 mass1 *= NewAxis1_E/((NewAxis1_E+(*Object[i]).energy())*(NewAxis1_E+(*Object[i]).energy()));
00296
00297 mass2 *= NewAxis2_E/((NewAxis2_E+(*Object[i]).energy())*(NewAxis2_E+(*Object[i]).energy()));
00298
00299 }
00300
00301 if(mass1 < mass2) {
00302 if (Object_Group[i] != 1){
00303 I_Move = true;
00304 }
00305 Object_Group[i] = 1;
00306
00307 Sum1_Px += (*Object[i]).px();
00308 Sum1_Py += (*Object[i]).py();
00309 Sum1_Pz += (*Object[i]).pz();
00310 Sum1_P += (*Object[i]).p();
00311 Sum1_E += (*Object[i]).energy();
00312 } else {
00313 if (Object_Group[i] != 2){
00314 I_Move = true;
00315 }
00316 Object_Group[i] = 2;
00317 Sum2_Px += (*Object[i]).px();
00318 Sum2_Py += (*Object[i]).py();
00319 Sum2_Pz += (*Object[i]).pz();
00320 Sum2_P += (*Object[i]).p();
00321 Sum2_E += (*Object[i]).energy();
00322
00323 }
00324
00325
00326 } else {
00327
00328 if (Object_Group[i] == 1){
00329 Sum1_Px += (*Object[i]).px();
00330 Sum1_Py += (*Object[i]).py();
00331 Sum1_Pz += (*Object[i]).pz();
00332 Sum1_P += (*Object[i]).p();
00333 Sum1_E += (*Object[i]).energy();
00334 }
00335 if (Object_Group[i] == 2){
00336 Sum2_Px += (*Object[i]).px();
00337 Sum2_Py += (*Object[i]).py();
00338 Sum2_Pz += (*Object[i]).pz();
00339 Sum2_P += (*Object[i]).p();
00340 Sum2_E += (*Object[i]).energy();
00341 }
00342
00343
00344
00345 }
00346
00347
00348 }
00349 }
00350 } else {
00351 throw cms::Exception("Configuration")
00352 << "Please give a valid hemisphere association method!";
00353 }
00354
00355
00356
00357 Axis1[3] = sqrt(Sum1_Px*Sum1_Px + Sum1_Py*Sum1_Py + Sum1_Pz*Sum1_Pz);
00358 if (Axis1[3]<0.0001) {
00359 edm::LogWarning("HemisphereAlgo") << "ZERO objects in group 1! ";
00360 } else {
00361 Axis1[0] = Sum1_Px / Axis1[3];
00362 Axis1[1] = Sum1_Py / Axis1[3];
00363 Axis1[2] = Sum1_Pz / Axis1[3];
00364 Axis1[4] = Sum1_E;
00365 }
00366
00367
00368
00369 Axis2[3] = sqrt(Sum2_Px*Sum2_Px + Sum2_Py*Sum2_Py + Sum2_Pz*Sum2_Pz);
00370 if (Axis2[3]<0.0001) {
00371 edm::LogWarning("HemisphereAlgo") << "ZERO objects in group 2!";
00372 } else {
00373 Axis2[0] = Sum2_Px / Axis2[3];
00374 Axis2[1] = Sum2_Py / Axis2[3];
00375 Axis2[2] = Sum2_Pz / Axis2[3];
00376 Axis2[4] = Sum2_E;
00377 }
00378
00379 LogDebug("HemisphereAlgo") << " Grouping = ";
00380 for (int i=0;i<vsize;i++){
00381 LogTrace("HemisphereAlgo") << " " << Object_Group[i];
00382 }
00383 LogTrace("HemisphereAlgo") << std::endl;
00384
00385 }
00386
00387
00388 status = 1;
00389 return 1;
00390 }
00391
00392
00393