cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
N/A
Central exclusive processes event generator
Toggle main menu visibility
Main Page
Related Pages
Packages
Package List
Package Members
All
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
x
y
Functions
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
x
y
Variables
Typedefs
Enumerations
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
a
b
c
d
e
f
g
h
i
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Variables
a
b
c
d
e
f
g
h
i
k
l
m
n
p
q
r
s
t
u
v
w
x
y
z
Typedefs
Enumerations
Enumerator
b
c
d
e
g
h
i
l
m
n
p
r
t
u
w
x
y
z
Related Symbols
d
g
h
o
u
v
Files
File List
File Members
All
_
b
c
d
h
o
p
r
s
Functions
Variables
Macros
_
b
c
d
p
r
s
▼
CepGen
Reference manual
Bibliography
►
Packages
►
Classes
▼
Files
▼
File List
▼
include
▼
CepGen
►
Cards
►
Core
►
Event
►
EventFilter
►
FormFactors
▼
Integration
►
FunctionalIntegrand.h
►
FunctionIntegrand.h
►
GridParameters.h
►
GSLIntegrator.h
►
Integrand.h
►
Integrator.h
►
ProcessIntegrand.h
►
Modules
►
PartonFluxes
►
Physics
►
Process
►
StructureFunctions
►
Utils
►
Generator.h
►
Version.h
►
CepGenHepMC2
►
CepGenHepMC3
►
CepGenHerwig6
►
CepGenMadGraph
►
CepGenPythia6
►
CepGenPythia8
►
CepGenPython
►
CepGenRoot
►
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
Loading...
Searching...
No Matches
GridParameters.h
Go to the documentation of this file.
1
/*
2
* CepGen: a central exclusive processes event generator
3
* Copyright (C) 2013-2025 Laurent Forthomme
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
* 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
#ifndef CepGen_Integration_GridParameters_h
20
#define CepGen_Integration_GridParameters_h
21
22
#include <vector>
23
24
namespace
cepgen::utils
{
25
class
RandomGenerator;
26
}
// namespace cepgen::utils
27
28
namespace
cepgen
{
30
class
GridParameters
{
31
public
:
33
explicit
GridParameters
(
size_t
mbin,
size_t
num_dimensions);
34
35
using
coord_t
= std::vector<unsigned short>;
36
37
void
dump
()
const
;
38
39
inline
size_t
size
()
const
{
return
coordinates_.size(); }
41
inline
const
coord_t
&
n
(
size_t
coord)
const
{
return
coordinates_.at(coord); }
42
inline
float
globalMax
()
const
{
return
f_max_global_; }
43
45
inline
float
maxValue
(
size_t
coord)
const
{
return
f_max_.at(coord); }
46
inline
double
maxValueDiff
()
const
{
return
f_max_diff_; }
47
inline
double
maxHistValue
()
const
{
return
f_max_old_; }
48
49
void
setValue
(
size_t
,
float
);
51
void
shoot
(
utils::RandomGenerator
& random_generator,
size_t
coordinate, std::vector<double>& out)
const
;
53
inline
size_t
numPoints
(
size_t
coordinate)
const
{
return
num_points_.at(coordinate); }
55
inline
void
increment
(
size_t
coordinate) { num_points_.at(coordinate)++; }
56
57
inline
bool
prepared
()
const
{
return
gen_prepared_; }
58
inline
void
setPrepared
(
bool
prepared
=
true
) { gen_prepared_ =
prepared
; }
59
61
inline
float
correctionValue
()
const
{
return
correction_; }
63
inline
void
setCorrectionValue
(
float
correction) { correction_ = correction; }
64
bool
correct
(
size_t
);
65
66
void
rescale
(
size_t
,
float
);
67
void
initCorrectionCycle
(
size_t
,
float
);
68
69
private
:
70
void
generateCoordinates(
coord_t
&,
size_t
)
const
;
71
72
const
size_t
mbin_;
73
const
double
inv_mbin_;
74
size_t
num_dimensions_{0};
75
bool
gen_prepared_{
false
};
76
float
correction_{0.};
77
float
correction2_{0.};
78
std::vector<coord_t> coordinates_;
79
std::vector<size_t> num_points_;
80
std::vector<float> f_max_;
81
float
f_max_global_{0.};
82
float
f_max2_{0.};
83
float
f_max_diff_{0.};
84
float
f_max_old_{0.};
85
};
86
}
// namespace cepgen
87
88
#endif
58
inline
void
setPrepared
(
bool
prepared
=
true
) { gen_prepared_ =
prepared
; } {
…
}
30
class
GridParameters
{
…
};
cepgen::GridParameters
A parameters placeholder for the grid integration helper.
Definition
GridParameters.h:30
cepgen::GridParameters::rescale
void rescale(size_t, float)
cepgen::GridParameters::setValue
void setValue(size_t, float)
Set the function value for a given grid coordinate Shoot a phase space point for a grid coordinate.
cepgen::GridParameters::maxValueDiff
double maxValueDiff() const
Definition
GridParameters.h:46
cepgen::GridParameters::size
size_t size() const
Grid multiplicity Number of times a phase space point has been randomly selected.
Definition
GridParameters.h:39
cepgen::GridParameters::GridParameters
GridParameters(size_t mbin, size_t num_dimensions)
Build a generation grid for a ndim-dimensional phase space.
cepgen::GridParameters::increment
void increment(size_t coordinate)
Specify a new trial has been attempted for bin.
Definition
GridParameters.h:55
cepgen::GridParameters::initCorrectionCycle
void initCorrectionCycle(size_t, float)
cepgen::GridParameters::setPrepared
void setPrepared(bool prepared=true)
Mark the grid as prepared.
Definition
GridParameters.h:58
cepgen::GridParameters::maxHistValue
double maxHistValue() const
Definition
GridParameters.h:47
cepgen::GridParameters::prepared
bool prepared() const
Has the grid been prepared?
Definition
GridParameters.h:57
cepgen::GridParameters::setCorrectionValue
void setCorrectionValue(float correction)
Set the correction to apply on the next phase space point generation.
Definition
GridParameters.h:63
cepgen::GridParameters::correct
bool correct(size_t)
Apply the correction requested at the previous generation.
cepgen::GridParameters::maxValue
float maxValue(size_t coord) const
Maximal function value for a given grid coordinate.
Definition
GridParameters.h:45
cepgen::GridParameters::n
const coord_t & n(size_t coord) const
Definition
GridParameters.h:41
cepgen::GridParameters::correctionValue
float correctionValue() const
Correction to apply on the next phase space point generation.
Definition
GridParameters.h:61
cepgen::GridParameters::coord_t
std::vector< unsigned short > coord_t
Coordinates definition.
Definition
GridParameters.h:35
cepgen::GridParameters::numPoints
size_t numPoints(size_t coordinate) const
Number of points already shot for a given grid coordinate.
Definition
GridParameters.h:53
cepgen::GridParameters::shoot
void shoot(utils::RandomGenerator &random_generator, size_t coordinate, std::vector< double > &out) const
cepgen::GridParameters::dump
void dump() const
Dump the grid coordinates.
cepgen::GridParameters::globalMax
float globalMax() const
Global function maximum.
Definition
GridParameters.h:42
cepgen::utils::RandomGenerator
A random number generator.
Definition
RandomGenerator.h:28
cepgen::utils
Collection of utilities.
Definition
RunParameters.h:33
cepgen
Common namespace for this Monte Carlo generator.
Definition
Handler.h:26
include
CepGen
Integration
GridParameters.h
Generated on Tue Apr 22 2025 for CepGen by
1.10.0