cepgen
is hosted by
Hepforge
,
IPPP Durham
CepGen
1.2.5
Central exclusive processes event generator
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-2023 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 <cstddef>
23
#include <vector>
24
25
namespace
cepgen
{
26
class
Integrator;
28
class
GridParameters
{
29
public
:
31
explicit
GridParameters
(
size_t
mbin,
size_t
ndim);
32
33
typedef
std::vector<unsigned short>
coord_t
;
34
35
void
dump
()
const
;
36
37
inline
size_t
size
()
const
{
return
coords_.size(); }
39
inline
const
coord_t
&
n
(
size_t
coord)
const
{
return
coords_.at(coord); }
40
inline
float
globalMax
()
const
{
return
f_max_global_; }
41
43
inline
float
maxValue
(
size_t
coord)
const
{
return
f_max_.at(coord); }
44
inline
double
maxValueDiff
()
const
{
return
f_max_diff_; }
45
inline
double
maxHistValue
()
const
{
return
f_max_old_; }
46
47
void
setValue
(
size_t
,
float
);
49
void
shoot
(
const
Integrator
* integ,
size_t
coord, std::vector<double>& out)
const
;
51
inline
size_t
numPoints
(
size_t
coord)
const
{
return
num_points_.at(coord); }
53
inline
void
increment
(
size_t
coord) { num_points_.at(coord)++; }
54
55
inline
bool
prepared
()
const
{
return
gen_prepared_; }
56
inline
void
setPrepared
(
bool
prepared
=
true
) { gen_prepared_ =
prepared
; }
57
59
inline
float
correctionValue
()
const
{
return
correc_; }
61
inline
void
setCorrectionValue
(
float
correc) { correc_ = correc; }
62
bool
correct
(
size_t
);
63
64
void
rescale
(
size_t
,
float
);
65
void
initCorrectionCycle
(
size_t
,
float
);
66
67
private
:
68
void
generateCoordinates(
coord_t
&,
size_t
)
const
;
69
70
const
size_t
mbin_;
71
const
double
inv_mbin_;
72
size_t
ndim_{0};
73
bool
gen_prepared_{
false
};
74
float
correc_{0.};
75
float
correc2_{0.};
76
std::vector<coord_t> coords_;
77
std::vector<size_t> num_points_;
78
std::vector<float> f_max_;
79
float
f_max_global_{0.};
80
float
f_max2_{0.};
81
float
f_max_diff_{0.};
82
float
f_max_old_{0.};
83
};
84
}
// namespace cepgen
85
86
#endif
cepgen::GridParameters
A parameters placeholder for the grid integration helper.
Definition
GridParameters.h:28
cepgen::GridParameters::rescale
void rescale(size_t, float)
Definition
GridParameters.cpp:86
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.
Definition
GridParameters.cpp:39
cepgen::GridParameters::maxValueDiff
double maxValueDiff() const
Definition
GridParameters.h:44
cepgen::GridParameters::size
size_t size() const
Grid multiplicity Number of times a phase space point has been randomly selected.
Definition
GridParameters.h:37
cepgen::GridParameters::increment
void increment(size_t coord)
Specify a new trial has been attempted for bin.
Definition
GridParameters.h:53
cepgen::GridParameters::initCorrectionCycle
void initCorrectionCycle(size_t, float)
Definition
GridParameters.cpp:94
cepgen::GridParameters::setPrepared
void setPrepared(bool prepared=true)
Mark the grid as prepared.
Definition
GridParameters.h:56
cepgen::GridParameters::maxHistValue
double maxHistValue() const
Definition
GridParameters.h:45
cepgen::GridParameters::prepared
bool prepared() const
Has the grid been prepared.
Definition
GridParameters.h:55
cepgen::GridParameters::correct
bool correct(size_t)
Apply the correction requested at the previous generation.
Definition
GridParameters.cpp:71
cepgen::GridParameters::maxValue
float maxValue(size_t coord) const
Maximal function value for a given grid coordinate.
Definition
GridParameters.h:43
cepgen::GridParameters::numPoints
size_t numPoints(size_t coord) const
Number of points already shot for a given grid coordinate.
Definition
GridParameters.h:51
cepgen::GridParameters::n
const coord_t & n(size_t coord) const
Definition
GridParameters.h:39
cepgen::GridParameters::correctionValue
float correctionValue() const
Correction to apply on the next phase space point generation.
Definition
GridParameters.h:59
cepgen::GridParameters::coord_t
std::vector< unsigned short > coord_t
Coordinates definition.
Definition
GridParameters.h:33
cepgen::GridParameters::shoot
void shoot(const Integrator *integ, size_t coord, std::vector< double > &out) const
Definition
GridParameters.cpp:45
cepgen::GridParameters::setCorrectionValue
void setCorrectionValue(float correc)
Set the correction to apply on the next phase space point generation.
Definition
GridParameters.h:61
cepgen::GridParameters::dump
void dump() const
Dump the grid coordinates.
Definition
GridParameters.cpp:52
cepgen::GridParameters::globalMax
float globalMax() const
Global function maximum.
Definition
GridParameters.h:40
cepgen::Integrator
Monte-Carlo integration algorithm.
Definition
Integrator.h:28
cepgen
Common namespace for this Monte Carlo generator.
Definition
CommandLineHandler.cpp:36
CepGen
Integration
GridParameters.h
Generated on Mon Jul 29 2024 for CepGen by
1.9.7