CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
MixBoostEvtVtxGenerator Class Reference
Inheritance diagram for MixBoostEvtVtxGenerator:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void Alpha (double m=0)
 angle between crossing plane and horizontal plane More...
 
void Beta (double m=0)
 
double BetaFunction (double z, double z0)
 beta function More...
 
void betastar (double m=0)
 set beta_star More...
 
void emittance (double m=0)
 emittance (no the normalized) More...
 
CLHEP::HepRandomEngine & getEngine ()
 
virtual TMatrixD * GetInvLorentzBoost ()
 
virtual HepMC::FourVector * getRecVertex (edm::Event &)
 
virtual HepMC::FourVector * getVertex (edm::Event &)
 
 MixBoostEvtVtxGenerator (const edm::ParameterSet &p)
 
virtual HepMC::FourVector * newVertex ()
 return a new event vertex More...
 
void Phi (double m=0)
 set half crossing angle More...
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 
void sigmaZ (double s=1.0)
 set resolution in Z in cm More...
 
void X0 (double m=0)
 set mean in X in cm More...
 
void Y0 (double m=0)
 set mean in Y in cm More...
 
void Z0 (double m=0)
 set mean in Z in cm More...
 
virtual ~MixBoostEvtVtxGenerator ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

 MixBoostEvtVtxGenerator (const MixBoostEvtVtxGenerator &p)
 
MixBoostEvtVtxGeneratoroperator= (const MixBoostEvtVtxGenerator &rhs)
 

Private Attributes

double alpha_
 
double beta_
 
TMatrixD * boost_
 
double falpha
 
double fbetastar
 
double femittance
 
CLHEP::HepRandomEngine * fEngine
 
CLHEP::RandGaussQ * fRandom
 
double fSigmaZ
 
double fTimeOffset
 
HepMC::FourVector * fVertex
 
double fX0
 
double fY0
 
double fZ0
 
edm::InputTag hiLabel
 
double phi_
 
edm::InputTag signalLabel
 
edm::InputTag sourceLabel
 
bool useRecVertex
 
std::vector< double > vtxOffset
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 52 of file MixBoostEvtVtxGenerator.cc.

Constructor & Destructor Documentation

MixBoostEvtVtxGenerator::MixBoostEvtVtxGenerator ( const edm::ParameterSet p)

Definition at line 122 of file MixBoostEvtVtxGenerator.cc.

References edm::hlt::Exception, edm::ParameterSet::exists(), fEngine, edm::RandomNumberGenerator::getEngine(), edm::ParameterSet::getParameter(), edm::Service< T >::isAvailable(), and vtxOffset.

122  :
123  fVertex(0), boost_(0), fTimeOffset(0), fEngine(0),
124  signalLabel(pset.getParameter<edm::InputTag>("signalLabel")),
125  hiLabel(pset.getParameter<edm::InputTag>("heavyIonLabel")),
126  useRecVertex(pset.exists("useRecVertex")?pset.getParameter<bool>("useRecVertex"):false)
127 {
128 
129  vtxOffset.resize(3);
130  if(pset.exists("vtxOffset")) vtxOffset=pset.getParameter< std::vector<double> >("vtxOffset");
132 
133  if ( ! rng.isAvailable()) {
134 
135  throw cms::Exception("Configuration")
136  << "The BaseEvtVtxGenerator requires the RandomNumberGeneratorService\n"
137  "which is not present in the configuration file. You must add the service\n"
138  "in the configuration file or remove the modules that require it.";
139  }
140 
141  CLHEP::HepRandomEngine& engine = rng->getEngine();
142  fEngine = &engine;
143 
144  produces<bool>("matchedVertex");
145 
146 }
CLHEP::HepRandomEngine * fEngine
std::vector< double > vtxOffset
bool isAvailable() const
Definition: Service.h:47
virtual CLHEP::HepRandomEngine & getEngine() const =0
Use this to get the random number engine, this is the only function most users should call...
MixBoostEvtVtxGenerator::~MixBoostEvtVtxGenerator ( )
virtual

Definition at line 148 of file MixBoostEvtVtxGenerator.cc.

References boost_, fRandom, and fVertex.

149 {
150  delete fVertex ;
151  if (boost_ != 0 ) delete boost_;
152  delete fRandom;
153 }
MixBoostEvtVtxGenerator::MixBoostEvtVtxGenerator ( const MixBoostEvtVtxGenerator p)
private

Copy constructor

Member Function Documentation

void MixBoostEvtVtxGenerator::Alpha ( double  m = 0)
inline

angle between crossing plane and horizontal plane

Definition at line 78 of file MixBoostEvtVtxGenerator.cc.

References m.

void MixBoostEvtVtxGenerator::Beta ( double  m = 0)
inline

Definition at line 79 of file MixBoostEvtVtxGenerator.cc.

References m.

double MixBoostEvtVtxGenerator::BetaFunction ( double  z,
double  z0 
)

beta function

Definition at line 188 of file MixBoostEvtVtxGenerator.cc.

References fbetastar, femittance, and mathSSE::sqrt().

Referenced by newVertex().

189 {
190  return sqrt(femittance*(fbetastar+(((z-z0)*(z-z0))/fbetastar)));
191 
192 }
float float float z
T sqrt(T t)
Definition: SSEVec.h:48
void MixBoostEvtVtxGenerator::betastar ( double  m = 0)
inline

set beta_star

Definition at line 82 of file MixBoostEvtVtxGenerator.cc.

References m.

void MixBoostEvtVtxGenerator::emittance ( double  m = 0)
inline

emittance (no the normalized)

Definition at line 84 of file MixBoostEvtVtxGenerator.cc.

References m.

CLHEP::HepRandomEngine & MixBoostEvtVtxGenerator::getEngine ( )

Definition at line 155 of file MixBoostEvtVtxGenerator.cc.

References fEngine.

155  {
156  return *fEngine;
157 }
CLHEP::HepRandomEngine * fEngine
TMatrixD * MixBoostEvtVtxGenerator::GetInvLorentzBoost ( )
virtual

Definition at line 207 of file MixBoostEvtVtxGenerator.cc.

References alpha_, beta_, boost_, funct::cos(), phi_, funct::sin(), mathSSE::sqrt(), and funct::tan().

207  {
208 
209  //alpha_ = 0;
210  //phi_ = 142.e-6;
211 // if (boost_ != 0 ) return boost_;
212 
213  //boost_.ResizeTo(4,4);
214  //boost_ = new TMatrixD(4,4);
215  TMatrixD tmpboost(4,4);
216  TMatrixD tmpboostZ(4,4);
217  TMatrixD tmpboostXYZ(4,4);
218 
219  //if ( (alpha_ == 0) && (phi_==0) ) { boost_->Zero(); return boost_; }
220 
221  // Lorentz boost to frame where the collision is head-on
222  // phi is the half crossing angle in the plane ZS
223  // alpha is the angle to the S axis from the X axis in the XY plane
224 
225  tmpboost(0,0) = 1./cos(phi_);
226  tmpboost(0,1) = - cos(alpha_)*sin(phi_);
227  tmpboost(0,2) = - tan(phi_)*sin(phi_);
228  tmpboost(0,3) = - sin(alpha_)*sin(phi_);
229  tmpboost(1,0) = - cos(alpha_)*tan(phi_);
230  tmpboost(1,1) = 1.;
231  tmpboost(1,2) = cos(alpha_)*tan(phi_);
232  tmpboost(1,3) = 0.;
233  tmpboost(2,0) = 0.;
234  tmpboost(2,1) = - cos(alpha_)*sin(phi_);
235  tmpboost(2,2) = cos(phi_);
236  tmpboost(2,3) = - sin(alpha_)*sin(phi_);
237  tmpboost(3,0) = - sin(alpha_)*tan(phi_);
238  tmpboost(3,1) = 0.;
239  tmpboost(3,2) = sin(alpha_)*tan(phi_);
240  tmpboost(3,3) = 1.;
241  //cout<<"beta "<<beta_;
242  double gama=1.0/sqrt(1-beta_*beta_);
243  tmpboostZ(0,0)=gama;
244  tmpboostZ(0,1)=0.;
245  tmpboostZ(0,2)=-1.0*beta_*gama;
246  tmpboostZ(0,3)=0.;
247  tmpboostZ(1,0)=0.;
248  tmpboostZ(1,1) = 1.;
249  tmpboostZ(1,2)=0.;
250  tmpboostZ(1,3)=0.;
251  tmpboostZ(2,0)=-1.0*beta_*gama;
252  tmpboostZ(2,1) = 0.;
253  tmpboostZ(2,2)=gama;
254  tmpboostZ(2,3) = 0.;
255  tmpboostZ(3,0)=0.;
256  tmpboostZ(3,1)=0.;
257  tmpboostZ(3,2)=0.;
258  tmpboostZ(3,3) = 1.;
259 
260  tmpboostXYZ=tmpboost*tmpboostZ;
261  tmpboost.Invert();
262 
263 
264 
265  boost_ = new TMatrixD(tmpboostXYZ);
266  boost_->Print();
267 
268  return boost_;
269 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
HepMC::FourVector * MixBoostEvtVtxGenerator::getRecVertex ( edm::Event evt)
virtual

Definition at line 304 of file MixBoostEvtVtxGenerator.cc.

References fVertex, edm::Event::getByLabel(), hiLabel, LaserDQM_cfg::input, and vtxOffset.

Referenced by produce().

304  {
305 
307  evt.getByLabel(hiLabel,input);
308 
309  double aX,aY,aZ;
310 
311  aX = input->begin()->position().x() + vtxOffset[0];
312  aY = input->begin()->position().y() + vtxOffset[1];
313  aZ = input->begin()->position().z() + vtxOffset[2];
314 
315  if(!fVertex) fVertex = new HepMC::FourVector();
316  fVertex->set(10.0*aX,10.0*aY,10.0*aZ,0.0); // HepMC positions in mm (RECO in cm)
317 
318  return fVertex;
319 
320 }
std::vector< double > vtxOffset
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
HepMC::FourVector * MixBoostEvtVtxGenerator::getVertex ( edm::Event evt)
virtual

Definition at line 271 of file MixBoostEvtVtxGenerator.cc.

References gather_cfg::cout, fVertex, edm::Event::getByLabel(), hiLabel, and LaserDQM_cfg::input.

Referenced by produce().

271  {
272 
274  evt.getByLabel(hiLabel,input);
275 
276  const HepMC::GenEvent* inev = input->GetEvent();
277  HepMC::GenVertex* genvtx = inev->signal_process_vertex();
278  if(!genvtx){
279  cout<<"No Signal Process Vertex!"<<endl;
280  HepMC::GenEvent::particle_const_iterator pt=inev->particles_begin();
281  HepMC::GenEvent::particle_const_iterator ptend=inev->particles_end();
282  while(!genvtx || ( genvtx->particles_in_size() == 1 && pt != ptend ) ){
283  if(!genvtx) cout<<"No Gen Vertex!"<<endl;
284  if(pt == ptend) cout<<"End reached!"<<endl;
285  genvtx = (*pt)->production_vertex();
286  ++pt;
287  }
288  }
289  double aX,aY,aZ,aT;
290 
291  aX = genvtx->position().x();
292  aY = genvtx->position().y();
293  aZ = genvtx->position().z();
294  aT = genvtx->position().t();
295 
296  if(!fVertex) fVertex = new HepMC::FourVector();
297  fVertex->set(aX,aY,aZ,aT);
298 
299  return fVertex;
300 
301 }
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
tuple cout
Definition: gather_cfg.py:121
HepMC::FourVector * MixBoostEvtVtxGenerator::newVertex ( )
virtual

return a new event vertex

Definition at line 161 of file MixBoostEvtVtxGenerator.cc.

References BetaFunction(), fRandom, fSigmaZ, fTimeOffset, fVertex, fX0, fY0, fZ0, mathSSE::sqrt(), X, and Gflash::Z.

161  {
162 
163 
164  double X,Y,Z;
165 
166  double tmp_sigz = fRandom->fire(0., fSigmaZ);
167  Z = tmp_sigz + fZ0;
168 
169  double tmp_sigx = BetaFunction(Z,fZ0);
170  // need sqrt(2) for beamspot width relative to single beam width
171  tmp_sigx /= sqrt(2.0);
172  X = fRandom->fire(0.,tmp_sigx) + fX0; // + Z*fdxdz ;
173 
174  double tmp_sigy = BetaFunction(Z,fZ0);
175  // need sqrt(2) for beamspot width relative to single beam width
176  tmp_sigy /= sqrt(2.0);
177  Y = fRandom->fire(0.,tmp_sigy) + fY0; // + Z*fdydz;
178 
179  double tmp_sigt = fRandom->fire(0., fSigmaZ);
180  double T = tmp_sigt + fTimeOffset;
181 
182  if ( fVertex == 0 ) fVertex = new HepMC::FourVector();
183  fVertex->set(X,Y,Z,T);
184 
185  return fVertex;
186 }
const double Z[kNumberCalorimeter]
#define X(str)
Definition: MuonsGrabber.cc:49
T sqrt(T t)
Definition: SSEVec.h:48
long double T
double BetaFunction(double z, double z0)
beta function
MixBoostEvtVtxGenerator& MixBoostEvtVtxGenerator::operator= ( const MixBoostEvtVtxGenerator rhs)
private

Copy assignment operator

void MixBoostEvtVtxGenerator::Phi ( double  m = 0)
inline

set half crossing angle

Definition at line 76 of file MixBoostEvtVtxGenerator.cc.

References m.

void MixBoostEvtVtxGenerator::produce ( edm::Event evt,
const edm::EventSetup  
)
virtual

Implements edm::EDProducer.

Definition at line 323 of file MixBoostEvtVtxGenerator.cc.

References edm::Event::getByLabel(), getRecVertex(), getVertex(), edm::Event::put(), hitfit::return, signalLabel, and useRecVertex.

324 {
325 
326 
327  Handle<HepMCProduct> HepMCEvt ;
328 
329  evt.getByLabel( signalLabel, HepMCEvt ) ;
330 
331  // generate new vertex & apply the shift
332  //
333 
334  HepMCEvt->applyVtxGen( useRecVertex ? getRecVertex(evt) : getVertex(evt) ) ;
335 
336  // HepMCEvt->boostToLab( GetInvLorentzBoost(), "vertex" );
337  // HepMCEvt->boostToLab( GetInvLorentzBoost(), "momentum" );
338 
339  // OK, create a (pseudo)product and put in into edm::Event
340  //
341  auto_ptr<bool> NewProduct(new bool(true)) ;
342  evt.put( NewProduct ,"matchedVertex") ;
343 
344  return ;
345 
346 }
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
virtual HepMC::FourVector * getRecVertex(edm::Event &)
virtual HepMC::FourVector * getVertex(edm::Event &)
void MixBoostEvtVtxGenerator::sigmaZ ( double  s = 1.0)

set resolution in Z in cm

Definition at line 195 of file MixBoostEvtVtxGenerator.cc.

References edm::hlt::Exception, fSigmaZ, and alignCSCRings::s.

196 {
197  if (s>=0 ) {
198  fSigmaZ=s;
199  }
200  else {
201  throw cms::Exception("LogicError")
202  << "Error in MixBoostEvtVtxGenerator::sigmaZ: "
203  << "Illegal resolution in Z (negative)";
204  }
205 }
void MixBoostEvtVtxGenerator::X0 ( double  m = 0)
inline

set mean in X in cm

Definition at line 69 of file MixBoostEvtVtxGenerator.cc.

References m.

void MixBoostEvtVtxGenerator::Y0 ( double  m = 0)
inline

set mean in Y in cm

Definition at line 71 of file MixBoostEvtVtxGenerator.cc.

References m.

void MixBoostEvtVtxGenerator::Z0 ( double  m = 0)
inline

set mean in Z in cm

Definition at line 73 of file MixBoostEvtVtxGenerator.cc.

References m.

Member Data Documentation

double MixBoostEvtVtxGenerator::alpha_
private

Definition at line 96 of file MixBoostEvtVtxGenerator.cc.

Referenced by GetInvLorentzBoost().

double MixBoostEvtVtxGenerator::beta_
private

Definition at line 98 of file MixBoostEvtVtxGenerator.cc.

Referenced by GetInvLorentzBoost().

TMatrixD* MixBoostEvtVtxGenerator::boost_
private

Definition at line 106 of file MixBoostEvtVtxGenerator.cc.

Referenced by GetInvLorentzBoost(), and ~MixBoostEvtVtxGenerator().

double MixBoostEvtVtxGenerator::falpha
private

Definition at line 103 of file MixBoostEvtVtxGenerator.cc.

double MixBoostEvtVtxGenerator::fbetastar
private

Definition at line 102 of file MixBoostEvtVtxGenerator.cc.

Referenced by BetaFunction().

double MixBoostEvtVtxGenerator::femittance
private

Definition at line 102 of file MixBoostEvtVtxGenerator.cc.

Referenced by BetaFunction().

CLHEP::HepRandomEngine* MixBoostEvtVtxGenerator::fEngine
private

Definition at line 109 of file MixBoostEvtVtxGenerator.cc.

Referenced by getEngine(), and MixBoostEvtVtxGenerator().

CLHEP::RandGaussQ* MixBoostEvtVtxGenerator::fRandom
private

Definition at line 112 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex(), and ~MixBoostEvtVtxGenerator().

double MixBoostEvtVtxGenerator::fSigmaZ
private

Definition at line 100 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex(), and sigmaZ().

double MixBoostEvtVtxGenerator::fTimeOffset
private

Definition at line 107 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex().

HepMC::FourVector* MixBoostEvtVtxGenerator::fVertex
private
double MixBoostEvtVtxGenerator::fX0
private

Definition at line 99 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex().

double MixBoostEvtVtxGenerator::fY0
private

Definition at line 99 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex().

double MixBoostEvtVtxGenerator::fZ0
private

Definition at line 99 of file MixBoostEvtVtxGenerator.cc.

Referenced by newVertex().

edm::InputTag MixBoostEvtVtxGenerator::hiLabel
private

Definition at line 115 of file MixBoostEvtVtxGenerator.cc.

Referenced by getRecVertex(), and getVertex().

double MixBoostEvtVtxGenerator::phi_
private

Definition at line 96 of file MixBoostEvtVtxGenerator.cc.

Referenced by GetInvLorentzBoost().

edm::InputTag MixBoostEvtVtxGenerator::signalLabel
private

Definition at line 114 of file MixBoostEvtVtxGenerator.cc.

Referenced by produce().

edm::InputTag MixBoostEvtVtxGenerator::sourceLabel
private

Definition at line 110 of file MixBoostEvtVtxGenerator.cc.

bool MixBoostEvtVtxGenerator::useRecVertex
private

Definition at line 116 of file MixBoostEvtVtxGenerator.cc.

Referenced by produce().

std::vector<double> MixBoostEvtVtxGenerator::vtxOffset
private

Definition at line 117 of file MixBoostEvtVtxGenerator.cc.

Referenced by getRecVertex(), and MixBoostEvtVtxGenerator().