forked from cepgen/cepgen
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Generator.cpp
219 lines (179 loc) · 8.29 KB
/
Generator.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
/*
* CepGen: a central exclusive processes event generator
* Copyright (C) 2013-2023 Laurent Forthomme
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*/
#include <chrono>
#include "CepGen/Core/Exception.h"
#include "CepGen/Core/GeneratorWorker.h"
#include "CepGen/EventFilter/EventExporter.h"
#include "CepGen/EventFilter/EventModifier.h"
#include "CepGen/Generator.h"
#include "CepGen/Integration/GridParameters.h"
#include "CepGen/Integration/Integrator.h"
#include "CepGen/Integration/ProcessIntegrand.h"
#include "CepGen/Modules/GeneratorWorkerFactory.h"
#include "CepGen/Modules/IntegratorFactory.h"
#include "CepGen/Parameters.h"
#include "CepGen/Process/Process.h"
#include "CepGen/Utils/String.h"
#include "CepGen/Utils/TimeKeeper.h"
namespace cepgen {
Generator::Generator(bool safe_mode) : parameters_(new Parameters) {
static bool init = false;
if (!init) {
cepgen::initialise(safe_mode);
init = true;
CG_DEBUG("Generator:init") << "Generator initialised";
}
//--- random number initialization
std::chrono::system_clock::time_point time = std::chrono::system_clock::now();
srandom(time.time_since_epoch().count());
}
Generator::Generator(Parameters* ip) : parameters_(ip) {}
Generator::~Generator() {
if (parameters_->timeKeeper())
CG_INFO("Generator:destructor") << parameters_->timeKeeper()->summary();
}
void Generator::clearRun() {
CG_DEBUG("Generator:clearRun") << "Run is set to be cleared.";
worker_ = GeneratorWorkerFactory::get().build(parameters_->generation().parameters().get<ParametersList>("worker"));
CG_DEBUG("Generator:clearRun") << "Initialised a generator worker with parameters: " << worker_->parameters()
<< ".";
// destroy and recreate the integrator instance
if (!integrator_)
resetIntegrator();
worker_->setRuntimeParameters(const_cast<const Parameters*>(parameters_.get()));
worker_->setIntegrator(integrator_.get());
xsect_ = Value{-1., -1.};
parameters_->prepareRun();
}
Parameters& Generator::parametersRef() { return *parameters_; }
void Generator::setParameters(Parameters* ip) { parameters_.reset(ip); }
double Generator::computePoint(const std::vector<double>& coord) {
if (!worker_)
clearRun();
if (!parameters_->hasProcess())
throw CG_FATAL("Generator:computePoint") << "Trying to compute a point with no process specified!";
const size_t ndim = worker_->integrand().process().ndim();
if (coord.size() != ndim)
throw CG_FATAL("Generator:computePoint") << "Invalid phase space dimension (ndim=" << ndim << ")!";
double res = worker_->integrand().eval(coord);
CG_DEBUG("Generator:computePoint") << "Result for x[" << ndim << "] = " << coord << ":\n\t" << res << ".";
return res;
}
void Generator::computeXsection(double& cross_section, double& err) {
const auto xsec = computeXsection();
cross_section = xsec;
err = xsec.uncertainty();
}
Value Generator::computeXsection() {
CG_INFO("Generator") << "Starting the computation of the process cross-section.";
integrate(); // run is cleared here
if (xsect_ < 1.e-2)
CG_INFO("Generator") << "Total cross section: " << xsect_ * 1.e3 << " fb.";
else if (xsect_ < 0.5e3)
CG_INFO("Generator") << "Total cross section: " << xsect_ << " pb.";
else if (xsect_ < 0.5e6)
CG_INFO("Generator") << "Total cross section: " << xsect_ * 1.e-3 << " nb.";
else if (xsect_ < 0.5e9)
CG_INFO("Generator") << "Total cross section: " << xsect_ * 1.e-6 << " µb.";
else
CG_INFO("Generator") << "Total cross section: " << xsect_ * 1.e-9 << " mb.";
return xsect_;
}
void Generator::resetIntegrator() {
CG_TICKER(parameters_->timeKeeper());
// create a spec-defined integrator in the current scope
setIntegrator(IntegratorFactory::get().build(parameters_->par_integrator));
}
void Generator::setIntegrator(std::unique_ptr<Integrator> integ) {
CG_TICKER(parameters_->timeKeeper());
// copy the integrator instance (or create it if unspecified) in the current scope
if (!integ)
integ = IntegratorFactory::get().build(parameters_->par_integrator);
integrator_ = std::move(integ);
CG_INFO("Generator:integrator") << "Generator will use a " << integrator_->name() << "-type integrator.";
}
void Generator::integrate() {
CG_TICKER(parameters_->timeKeeper());
clearRun();
if (!parameters_->hasProcess())
throw CG_FATAL("Generator:integrate") << "Trying to integrate while no process is specified!";
const size_t ndim = worker_->integrand().process().ndim();
if (ndim == 0)
throw CG_FATAL("Generator:integrate") << "Invalid phase space dimension. "
<< "At least one integration variable is required!";
CG_DEBUG("Generator:integrate") << "New integrator instance created for " << ndim << "-dimensional integration.";
if (!integrator_)
throw CG_FATAL("Generator:integrate") << "No integrator object was declared for the generator!";
xsect_ = integrator_->integrate(worker_->integrand());
CG_DEBUG("Generator:integrate") << "Computed cross section: (" << xsect_ << ") pb.";
// now that the cross section has been computed, feed it to the event modification algorithms...
for (auto& mod : parameters_->eventModifiersSequence())
mod->setCrossSection(xsect_);
// ...and to the event storage algorithms
for (auto& mod : parameters_->eventExportersSequence())
mod->setCrossSection(xsect_);
}
const Event& Generator::generateOneEvent(Event::callback callback) { return next(callback); }
void Generator::initialise() {
if (initialised_)
return;
if (!parameters_)
throw CG_FATAL("Generator:generate") << "No steering parameters specified!";
CG_TICKER(parameters_->timeKeeper());
// if no worker is found, launch a new integration/run preparation
if (!worker_)
integrate();
// prepare the run parameters for event generation
parameters_->initialise();
initialised_ = true;
}
const Event& Generator::next(Event::callback callback) {
if (!worker_ || !initialised_)
initialise();
size_t num_try = 0;
while (!worker_->next(callback)) {
if (num_try++ > 5)
throw CG_FATAL("Generator:next") << "Failed to generate the next event!";
}
return worker_->integrand().process().event();
}
void Generator::generate(size_t num_events, Event::callback callback) {
CG_TICKER(parameters_->timeKeeper());
initialise();
//--- if invalid argument, retrieve from runtime parameters
if (num_events < 1) {
if (parameters_->generation().targetLuminosity() > 0.) {
num_events = std::ceil(parameters_->generation().targetLuminosity() * xsect_);
CG_INFO("Generator") << "Target luminosity: " << parameters_->generation().targetLuminosity() << " pb-1.";
} else
num_events = parameters_->generation().maxGen();
}
CG_INFO("Generator") << utils::s("event", num_events, true) << " will be generated.";
const utils::Timer tmr;
//--- launch the event generation
worker_->generate(num_events, callback);
const double gen_time_s = tmr.elapsed();
const double rate_ms =
(parameters_->numGeneratedEvents() > 0) ? gen_time_s / parameters_->numGeneratedEvents() * 1.e3 : 0.;
const double equiv_lumi = parameters_->numGeneratedEvents() / crossSection();
CG_INFO("Generator") << utils::s("event", parameters_->numGeneratedEvents()) << " generated in " << gen_time_s
<< " s "
<< "(" << rate_ms << " ms/event).\n\t"
<< "Equivalent luminosity: " << utils::format("%g", equiv_lumi) << " pb^-1.";
}
} // namespace cepgen