Main Page
Namespaces
Classes
Package Documentation
IOMC
EventVertexGenerators
interface
BetafuncEvtVtxGenerator.h
Go to the documentation of this file.
1
#ifndef IOMC_BetafuncEvtVtxGenerator_H
2
#define IOMC_BetafuncEvtVtxGenerator_H
3
4
/*
5
________________________________________________________________________
6
7
BetafuncEvtVtxGenerator
8
9
Smear vertex according to the Beta function on the transverse plane
10
and a Gaussian on the z axis. It allows the beam to have a crossing
11
angle (dx/dz and dy/dz).
12
13
Based on GaussEvtVtxGenerator.h
14
implemented by Francisco Yumiceva (yumiceva@fnal.gov)
15
16
FERMILAB
17
2006
18
________________________________________________________________________
19
*/
20
21
#include "
IOMC/EventVertexGenerators/interface/BaseEvtVtxGenerator.h
"
22
#include "
FWCore/Framework/interface/ESWatcher.h
"
23
#include "
CondFormats/DataRecord/interface/SimBeamSpotObjectsRcd.h
"
24
25
namespace
CLHEP
{
26
class
HepRandomEngine;
27
}
28
29
class
BetafuncEvtVtxGenerator
:
public
BaseEvtVtxGenerator
30
{
31
public
:
32
BetafuncEvtVtxGenerator
(
const
edm::ParameterSet
&
p
);
33
virtual
~
BetafuncEvtVtxGenerator
();
34
35
virtual
void
beginRun(
const
edm::Run
& ,
const
edm::EventSetup
&)
override
;
36
virtual
void
beginLuminosityBlock(
edm::LuminosityBlock
const
&,
edm::EventSetup
const
&)
override
;
37
39
//virtual CLHEP::Hep3Vector * newVertex();
40
virtual
HepMC::FourVector* newVertex(CLHEP::HepRandomEngine*)
override
;
41
42
virtual
TMatrixD* GetInvLorentzBoost()
override
;
43
44
46
void
sigmaZ
(
double
s
=1.0);
47
49
void
X0
(
double
m
=0) { fX0=
m
; }
51
void
Y0
(
double
m
=0) { fY0=
m
; }
53
void
Z0
(
double
m
=0) { fZ0=
m
; }
54
56
void
Phi
(
double
m
=0) { phi_=
m
; }
58
void
Alpha
(
double
m
=0) { alpha_=
m
; }
59
61
void
betastar
(
double
m
=0) { fbetastar=
m
; }
63
void
emittance
(
double
m
=0) { femittance=
m
; }
64
66
double
BetaFunction(
double
z,
double
z0);
67
68
private
:
70
BetafuncEvtVtxGenerator
(
const
BetafuncEvtVtxGenerator
&p);
72
BetafuncEvtVtxGenerator
& operator = (
const
BetafuncEvtVtxGenerator
& rhs );
73
74
private
:
75
76
bool
readDB_
;
77
78
double
alpha_,
phi_
;
79
//TMatrixD boost_;
80
81
double
fX0, fY0,
fZ0
;
82
double
fSigmaZ
;
83
//double fdxdz, fdydz;
84
double
fbetastar,
femittance
;
85
// double falpha;
86
double
fTimeOffset
;
87
88
void
update
(
const
edm::EventSetup
& iEventSetup);
89
edm::ESWatcher<SimBeamSpotObjectsRcd>
parameterWatcher_
;
90
};
91
92
#endif
BetafuncEvtVtxGenerator::Phi
void Phi(double m=0)
set half crossing angle
Definition:
BetafuncEvtVtxGenerator.h:56
BetafuncEvtVtxGenerator::X0
void X0(double m=0)
set mean in X in cm
Definition:
BetafuncEvtVtxGenerator.h:49
ESWatcher.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition:
AlCaHLTBitMon_ParallelJobs.py:152
edm::LuminosityBlock
Definition:
LuminosityBlock.h:45
alignCSCRings.s
s
Definition:
alignCSCRings.py:91
SimBeamSpotObjectsRcd.h
BetafuncEvtVtxGenerator::fZ0
double fZ0
Definition:
BetafuncEvtVtxGenerator.h:81
BaseEvtVtxGenerator.h
fftjetvertexadder_cfi.sigmaZ
sigmaZ
Definition:
fftjetvertexadder_cfi.py:32
BetafuncEvtVtxGenerator::fTimeOffset
double fTimeOffset
Definition:
BetafuncEvtVtxGenerator.h:86
BetafuncEvtVtxGenerator::femittance
double femittance
Definition:
BetafuncEvtVtxGenerator.h:84
CLHEP
Definition:
CocoaGlobals.h:27
BetafuncEvtVtxGenerator
Definition:
BetafuncEvtVtxGenerator.h:29
BetafuncEvtVtxGenerator::fSigmaZ
double fSigmaZ
Definition:
BetafuncEvtVtxGenerator.h:82
BetafuncEvtVtxGenerator::betastar
void betastar(double m=0)
set beta_star
Definition:
BetafuncEvtVtxGenerator.h:61
BetafuncEvtVtxGenerator::Y0
void Y0(double m=0)
set mean in Y in cm
Definition:
BetafuncEvtVtxGenerator.h:51
edm::EventSetup
Definition:
EventSetup.h:45
edm::ESWatcher< SimBeamSpotObjectsRcd >
BetafuncEvtVtxGenerator::readDB_
bool readDB_
Definition:
BetafuncEvtVtxGenerator.h:76
BetafuncEvtVtxGenerator::parameterWatcher_
edm::ESWatcher< SimBeamSpotObjectsRcd > parameterWatcher_
Definition:
BetafuncEvtVtxGenerator.h:89
BaseEvtVtxGenerator
Definition:
BaseEvtVtxGenerator.h:26
BetafuncEvtVtxGenerator::Z0
void Z0(double m=0)
set mean in Z in cm
Definition:
BetafuncEvtVtxGenerator.h:53
BetafuncEvtVtxGenerator::Alpha
void Alpha(double m=0)
angle between crossing plane and horizontal plane
Definition:
BetafuncEvtVtxGenerator.h:58
funct::m
m
Definition:
Factorize.h:54
update
#define update(a, b)
Definition:
TrackClassifier.cc:10
BetafuncEvtVtxGenerator::phi_
double phi_
Definition:
BetafuncEvtVtxGenerator.h:78
edm::ParameterSet
Definition:
ParameterSet.h:36
BetafuncEvtVtxGenerator::emittance
void emittance(double m=0)
emittance (no the normalized)
Definition:
BetafuncEvtVtxGenerator.h:63
edm::Run
Definition:
Run.h:42
Generated for CMSSW Reference Manual by
1.8.11