Fork me on GitHub

source: git/modules/Weighter.cc@ e5ea42e

Last change on this file since e5ea42e was cab38f6, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove svn tags and fix formatting

  • Property mode set to 100644
File size: 4.7 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19
20/** \class Weighter
21 *
22 * Apply a weight depending on PDG code.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/Weighter.h"
29
30#include "classes/DelphesClasses.h"
31#include "classes/DelphesFactory.h"
32#include "classes/DelphesFormula.h"
33
34#include "ExRootAnalysis/ExRootResult.h"
35#include "ExRootAnalysis/ExRootFilter.h"
36#include "ExRootAnalysis/ExRootClassifier.h"
37
38#include "TMath.h"
39#include "TString.h"
40#include "TFormula.h"
41#include "TRandom3.h"
42#include "TObjArray.h"
43#include "TDatabasePDG.h"
44#include "TLorentzVector.h"
45
46#include <algorithm>
47#include <stdexcept>
48#include <iostream>
49#include <sstream>
50
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55bool Weighter::TIndexStruct::operator< (const Weighter::TIndexStruct &value) const
56{
57 Int_t i;
58
59 for(i = 0; i < 4; ++i)
60 {
61 if(codes[i] != value.codes[i]) return codes[i] < value.codes[i];
62 }
63
64 return false;
65}
66
67//------------------------------------------------------------------------------
68
69Weighter::Weighter() :
70 fItInputArray(0)
71{
72}
73
74//------------------------------------------------------------------------------
75
76Weighter::~Weighter()
77{
78}
79
80//------------------------------------------------------------------------------
81
82void Weighter::Init()
83{
84 ExRootConfParam param, paramCodes;
85 Int_t i, j, size, sizeCodes;
86 Int_t code;
87 TIndexStruct index;
88 Double_t weight;
89
90 fWeightSet.clear();
91
92
93 // set default weight value
94 fWeightMap.clear();
95 memset(index.codes, 0, 4*sizeof(Int_t));
96 fWeightMap[index] = 1.0;
97
98 // read weights
99 param = GetParam("Weight");
100 size = param.GetSize();
101 for(i = 0; i < size/2; ++i)
102 {
103 paramCodes = param[i*2];
104 sizeCodes = paramCodes.GetSize();
105 weight = param[i*2 + 1].GetDouble();
106
107 if(sizeCodes < 1 || sizeCodes > 4)
108 {
109 throw runtime_error("only 1, 2, 3 or 4 PDG codes can be specified per weight");
110 }
111
112 memset(index.codes, 0, 4*sizeof(Int_t));
113
114 for(j = 0; j < sizeCodes; ++j)
115 {
116 code = paramCodes[j].GetInt();
117 index.codes[j] = code;
118 fWeightSet.insert(code);
119 }
120
121 sort(index.codes, index.codes + 4);
122
123 fWeightMap[index] = weight;
124 }
125
126 // import input array(s)
127
128 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
129 fItInputArray = fInputArray->MakeIterator();
130
131 // create output array(s)
132
133 fOutputArray = ExportArray(GetString("OutputArray", "weight"));
134}
135
136//------------------------------------------------------------------------------
137
138void Weighter::Finish()
139{
140 if(fItInputArray) delete fItInputArray;
141}
142
143//------------------------------------------------------------------------------
144
145void Weighter::Process()
146{
147 Candidate *candidate;
148 Int_t i;
149 TIndexStruct index;
150 Double_t weight;
151 set<Int_t>::iterator itCodeSet;
152 map<TIndexStruct, Double_t>::iterator itWeightMap;
153
154 DelphesFactory *factory = GetFactory();
155
156 // loop over all particles
157 fCodeSet.clear();
158 fItInputArray->Reset();
159 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
160 {
161 if(candidate->Status != 3) continue;
162
163 if(fWeightSet.find(candidate->PID) == fWeightSet.end()) continue;
164
165 fCodeSet.insert(candidate->PID);
166 }
167
168 // find default weight value
169 memset(index.codes, 0, 4*sizeof(Int_t));
170 itWeightMap = fWeightMap.find(index);
171 weight = itWeightMap->second;
172
173 if(fCodeSet.size() <= 4)
174 {
175 i = 0;
176 for(itCodeSet = fCodeSet.begin(); itCodeSet != fCodeSet.end(); ++itCodeSet)
177 {
178 index.codes[i] = *itCodeSet;
179 ++i;
180 }
181
182 sort(index.codes, index.codes + 4);
183
184 itWeightMap = fWeightMap.find(index);
185 if(itWeightMap != fWeightMap.end())
186 {
187 weight = itWeightMap->second;
188 }
189 }
190
191 candidate = factory->NewCandidate();
192 candidate->Momentum.SetPtEtaPhiE(weight, 0.0, 0.0, weight);
193 fOutputArray->Add(candidate);
194}
195
196//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.