Fork me on GitHub

source: git/modules/Weighter.cc@ d836f45

ImprovedOutputFile Timing dual_readout llp
Last change on this file since d836f45 was 1fa50c2, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix GPLv3 header

  • 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 * $Date$
25 * $Revision$
26 *
27 *
28 * \author P. Demin - UCL, Louvain-la-Neuve
29 *
30 */
31
32#include "modules/Weighter.h"
33
34#include "classes/DelphesClasses.h"
35#include "classes/DelphesFactory.h"
36#include "classes/DelphesFormula.h"
37
38#include "ExRootAnalysis/ExRootResult.h"
39#include "ExRootAnalysis/ExRootFilter.h"
40#include "ExRootAnalysis/ExRootClassifier.h"
41
42#include "TMath.h"
43#include "TString.h"
44#include "TFormula.h"
45#include "TRandom3.h"
46#include "TObjArray.h"
47#include "TDatabasePDG.h"
48#include "TLorentzVector.h"
49
50#include <algorithm>
51#include <stdexcept>
52#include <iostream>
53#include <sstream>
54
55using namespace std;
56
57//------------------------------------------------------------------------------
58
59bool Weighter::TIndexStruct::operator< (const Weighter::TIndexStruct &value) const
60{
61 Int_t i;
62
63 for(i = 0; i < 4; ++i)
64 {
65 if(codes[i] != value.codes[i]) return codes[i] < value.codes[i];
66 }
67
68 return false;
69}
70
71//------------------------------------------------------------------------------
72
73Weighter::Weighter() :
74 fItInputArray(0)
75{
76}
77
78//------------------------------------------------------------------------------
79
80Weighter::~Weighter()
81{
82}
83
84//------------------------------------------------------------------------------
85
86void Weighter::Init()
87{
88 ExRootConfParam param, paramCodes;
89 Int_t i, j, size, sizeCodes;
90 Int_t code;
91 TIndexStruct index;
92 Double_t weight;
93
94 fWeightSet.clear();
95
96
97 // set default weight value
98 fWeightMap.clear();
99 memset(index.codes, 0, 4*sizeof(Int_t));
100 fWeightMap[index] = 1.0;
101
102 // read weights
103 param = GetParam("Weight");
104 size = param.GetSize();
105 for(i = 0; i < size/2; ++i)
106 {
107 paramCodes = param[i*2];
108 sizeCodes = paramCodes.GetSize();
109 weight = param[i*2 + 1].GetDouble();
110
111 if(sizeCodes < 1 || sizeCodes > 4)
112 {
113 throw runtime_error("only 1, 2, 3 or 4 PDG codes can be specified per weight");
114 }
115
116 memset(index.codes, 0, 4*sizeof(Int_t));
117
118 for(j = 0; j < sizeCodes; ++j)
119 {
120 code = paramCodes[j].GetInt();
121 index.codes[j] = code;
122 fWeightSet.insert(code);
123 }
124
125 sort(index.codes, index.codes + 4);
126
127 fWeightMap[index] = weight;
128 }
129
130 // import input array(s)
131
132 fInputArray = ImportArray(GetString("InputArray", "Delphes/allParticles"));
133 fItInputArray = fInputArray->MakeIterator();
134
135 // create output array(s)
136
137 fOutputArray = ExportArray(GetString("OutputArray", "weight"));
138}
139
140//------------------------------------------------------------------------------
141
142void Weighter::Finish()
143{
144 if(fItInputArray) delete fItInputArray;
145}
146
147//------------------------------------------------------------------------------
148
149void Weighter::Process()
150{
151 Candidate *candidate;
152 Int_t i;
153 TIndexStruct index;
154 Double_t weight;
155 set<Int_t>::iterator itCodeSet;
156 map<TIndexStruct, Double_t>::iterator itWeightMap;
157
158 DelphesFactory *factory = GetFactory();
159
160 // loop over all particles
161 fCodeSet.clear();
162 fItInputArray->Reset();
163 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
164 {
165 if(candidate->Status != 3) continue;
166
167 if(fWeightSet.find(candidate->PID) == fWeightSet.end()) continue;
168
169 fCodeSet.insert(candidate->PID);
170 }
171
172 // find default weight value
173 memset(index.codes, 0, 4*sizeof(Int_t));
174 itWeightMap = fWeightMap.find(index);
175 weight = itWeightMap->second;
176
177 if(fCodeSet.size() <= 4)
178 {
179 i = 0;
180 for(itCodeSet = fCodeSet.begin(); itCodeSet != fCodeSet.end(); ++itCodeSet)
181 {
182 index.codes[i] = *itCodeSet;
183 ++i;
184 }
185
186 sort(index.codes, index.codes + 4);
187
188 itWeightMap = fWeightMap.find(index);
189 if(itWeightMap != fWeightMap.end())
190 {
191 weight = itWeightMap->second;
192 }
193 }
194
195 candidate = factory->NewCandidate();
196 candidate->Momentum.SetPtEtaPhiE(weight, 0.0, 0.0, weight);
197 fOutputArray->Add(candidate);
198}
199
200//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.