CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

BeamProfileVtxGenerator Class Reference

#include <BeamProfileVtxGenerator.h>

Inheritance diagram for BeamProfileVtxGenerator:
BaseEvtVtxGenerator edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

void beamPos (double m=0)
 set mean in Z in cm
 BeamProfileVtxGenerator (const edm::ParameterSet &p)
void eta (double m=0)
 set eta
virtual TMatrixD * GetInvLorentzBoost ()
void meanX (double m=0)
 set mean in X in cm
void meanY (double m=0)
 set mean in Y in cm
virtual HepMC::FourVector * newVertex ()
 return a new event vertex
void phi (double m=0)
 set phi in radian
void psi (double m=999)
 set psi in radian
void setType (bool m=true)
 set type
void sigmaX (double s=1.0)
 set resolution in X in cm
void sigmaY (double s=1.0)
 set resolution in Y in cm
virtual ~BeamProfileVtxGenerator ()

Private Member Functions

 BeamProfileVtxGenerator (const BeamProfileVtxGenerator &p)
BeamProfileVtxGeneratoroperator= (const BeamProfileVtxGenerator &rhs)

Private Attributes

std::vector< double > fdistn
double fEta
bool ffile
double fMeanX
double fMeanY
double fMeanZ
double fPhi
double fPsi
CLHEP::HepRandom * fRandom
double fSigmaX
double fSigmaY
double fTheta
double fTimeOffset
bool fType
int nBinx
int nBiny

Detailed Description

Definition at line 21 of file BeamProfileVtxGenerator.h.


Constructor & Destructor Documentation

BeamProfileVtxGenerator::BeamProfileVtxGenerator ( const edm::ParameterSet p)

Definition at line 20 of file BeamProfileVtxGenerator.cc.

References beamPos(), HTMLExport::elem(), eta(), fdistn, fEta, ffile, mergeVDriftHistosByStation::file, fMeanX, fMeanY, fMeanZ, fPhi, fPsi, fSigmaX, fSigmaY, fTheta, fTimeOffset, fType, edm::ParameterSet::getParameter(), i, recoMuon::in, prof2calltree::last, meanX(), meanY(), nBinx, nBiny, phi(), psi(), setType(), sigmaX(), and sigmaY().

                                                                          :
  BaseEvtVtxGenerator(p), fRandom(0) {
  
  meanX(p.getParameter<double>("BeamMeanX")*cm);
  meanY(p.getParameter<double>("BeamMeanY")*cm);
  beamPos(p.getParameter<double>("BeamPosition")*cm);
  sigmaX(p.getParameter<double>("BeamSigmaX")*cm);
  sigmaY(p.getParameter<double>("BeamSigmaY")*cm);
  double fMinEta = p.getParameter<double>("MinEta");
  double fMaxEta = p.getParameter<double>("MaxEta");
  double fMinPhi = p.getParameter<double>("MinPhi");
  double fMaxPhi = p.getParameter<double>("MaxPhi");
  eta(0.5*(fMinEta+fMaxEta));
  phi(0.5*(fMinPhi+fMaxPhi));
  psi(p.getParameter<double>("Psi"));
  nBinx = p.getParameter<int>("BinX");
  nBiny = p.getParameter<int>("BinY");
  ffile = p.getParameter<bool>("UseFile");
  fTimeOffset = p.getParameter<double>("TimeOffset")*ns*c_light;
  
  if (ffile) {
    std::string file = p.getParameter<std::string>("File");
    ifstream is(file.c_str(), std::ios::in);
    if (is) {
      double elem,sum=0;
      while (!is.eof()) {
        is >> elem;
        fdistn.push_back(elem);
        sum += elem;
      }
      if (((int)(fdistn.size())) >= nBinx*nBiny) {
        double last = 0;
        for (unsigned int i=0; i<fdistn.size(); i++) {
          fdistn[i] /= sum;
          fdistn[i] += last;
          last       = fdistn[i];
        }
        setType(false);
      } else {
        ffile = false;
      }
    } else {
      ffile = false;
    }
  } 
  if (!ffile) {
    setType(p.getParameter<bool>("GaussianProfile"));
  }

  edm::LogInfo("VertexGenerator") << "BeamProfileVtxGenerator: with "
                                  << "beam along eta = " << fEta 
                                  << " (Theta = " << fTheta/deg 
                                  << ") phi = " << fPhi/deg 
                                  << ") psi = " << fPsi/deg 
                                  << " centred at (" << fMeanX << ", " 
                                  << fMeanY << ", "  << fMeanZ << ") "
                                  << "and spread (" << fSigmaX << ", "
                                  << fSigmaY << ") of type Gaussian = "
                                  << fType << " use file " << ffile;
  if (ffile) 
    edm::LogInfo("VertexGenerator") << "With " << nBinx << " bins "
                                    << " along X and " << nBiny 
                                    << " bins along Y";
}
BeamProfileVtxGenerator::~BeamProfileVtxGenerator ( ) [virtual]

Definition at line 85 of file BeamProfileVtxGenerator.cc.

References fRandom.

                                                  {
  delete fRandom;
}
BeamProfileVtxGenerator::BeamProfileVtxGenerator ( const BeamProfileVtxGenerator p) [private]

Copy constructor


Member Function Documentation

void BeamProfileVtxGenerator::beamPos ( double  m = 0) [inline]

set mean in Z in cm

Definition at line 46 of file BeamProfileVtxGenerator.h.

References fMeanZ, and m.

Referenced by BeamProfileVtxGenerator().

{fMeanZ=m;}
void BeamProfileVtxGenerator::eta ( double  m = 0)

set eta

Definition at line 197 of file BeamProfileVtxGenerator.cc.

References create_public_lumi_plots::exp, fEta, fTheta, and alignCSCRings::s.

Referenced by BeamProfileVtxGenerator().

                                          { 
  fEta   = s; 
  fTheta = 2.0*atan(exp(-fEta));
}
virtual TMatrixD* BeamProfileVtxGenerator::GetInvLorentzBoost ( ) [inline, virtual]

Implements BaseEvtVtxGenerator.

Definition at line 31 of file BeamProfileVtxGenerator.h.

                                         {
          return 0;
  }
void BeamProfileVtxGenerator::meanX ( double  m = 0) [inline]

set mean in X in cm

Definition at line 42 of file BeamProfileVtxGenerator.h.

References fMeanX, and m.

Referenced by BeamProfileVtxGenerator().

{fMeanX=m;}
void BeamProfileVtxGenerator::meanY ( double  m = 0) [inline]

set mean in Y in cm

Definition at line 44 of file BeamProfileVtxGenerator.h.

References fMeanY, and m.

Referenced by BeamProfileVtxGenerator().

{fMeanY=m;}
HepMC::FourVector * BeamProfileVtxGenerator::newVertex ( ) [virtual]

return a new event vertex

Implements BaseEvtVtxGenerator.

Definition at line 91 of file BeamProfileVtxGenerator.cc.

References funct::cos(), fdistn, ffile, fMeanX, fMeanY, fMeanZ, fPhi, fPsi, fRandom, fSigmaX, fSigmaY, fTheta, fTimeOffset, fType, BaseEvtVtxGenerator::fVertex, i, LogDebug, M_PI, nBinx, nBiny, diffTwoXMLs::r1, and funct::sin().

                                                    {
  double aX, aY;
  if (ffile) {
    double r1 = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire();
    int ixy = 0, ix, iy;
    for (unsigned int i=0; i<fdistn.size(); i++) {
      if (r1 > fdistn[i]) ixy = i+1;
    }
    if (ixy >= nBinx*nBiny) {
      ix = nBinx-1; iy = nBiny-1;
    } else {
      ix = ixy%nBinx; iy = (ixy-ix)/nBinx;
    }
    aX = 0.5*(2*ix-nBinx+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaX + fMeanX ;
    aY = 0.5*(2*iy-nBiny+2*(dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire())*fSigmaY + fMeanY ;
  } else {
    if (fType) {
      aX = fSigmaX*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanX;
      aY = fSigmaY*(dynamic_cast<CLHEP::RandGaussQ*>(fRandom))->fire() +fMeanY;
    } else {
      aX = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaX,0.5*fSigmaX) + fMeanX ;
      aY = (dynamic_cast<CLHEP::RandFlat*>(fRandom))->fire(-0.5*fSigmaY,0.5*fSigmaY) + fMeanY;
    }
  }

  double xp, yp, zp ;
  if( 2.*M_PI < fabs(fPsi) )
  {
     xp = -aX*cos(fTheta)*cos(fPhi) +aY*sin(fPhi) +fMeanZ*sin(fTheta)*cos(fPhi);
     yp = -aX*cos(fTheta)*sin(fPhi) -aY*cos(fPhi) +fMeanZ*sin(fTheta)*sin(fPhi);
     zp =  aX*sin(fTheta)                         +fMeanZ*cos(fTheta);
  }
  else // for endcap testbeam, use 3rd angle psi
  {
     const HepGeom::Vector3D<double>  av ( aX, aY, fMeanZ ) ;

/*
     static const double kRadToDeg ( 180./M_PI ) ;
     std::cout<<"theta = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<fTheta*kRadToDeg<<", phi="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              << fPhi*kRadToDeg <<", PSI="<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<fPsi*kRadToDeg<<std::endl ;
*/
     const HepGeom::RotateZ3D R1 ( fPhi - M_PI ) ;
     const HepGeom::Point3D<double>  xUnit ( 0,1,0 ) ;
     const HepGeom::Point3D<double>  zUnit ( 0,0,1 ) ;
     const HepGeom::Transform3D RXRZ ( HepGeom::Rotate3D( -fTheta, R1*xUnit )*R1 ) ;
     const HepGeom::Transform3D TRF ( HepGeom::Rotate3D(fPsi,RXRZ*zUnit)*RXRZ ) ;
/*
     std::cout<<"\n\n$$$$$$$$$$$Transform="
              <<" thetaZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().thetaZ()*kRadToDeg
              <<", phiZ = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().phiZ()*kRadToDeg
              <<", thetaY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().thetaY()*kRadToDeg
              <<", phiY = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().phiY()*kRadToDeg
              <<", thetaX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().thetaX()*kRadToDeg
              <<", phiX = "<<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getRotation().phiX()*kRadToDeg
              <<std::setw(7) << std::setiosflags( std::ios::fixed ) << std::setprecision(5)
              <<TRF.getTranslation()<<std::endl ;
*/
     const HepGeom::Vector3D<double>  pv ( TRF*av ) ;

     xp = pv.x() ;
     yp = pv.y() ;
     zp = pv.z() ;
  }
  //if (fVertex == 0) fVertex = new CLHEP::Hep3Vector;
  //fVertex->set(xp, yp, zp);
  if (fVertex == 0 ) fVertex = new HepMC::FourVector() ;
  fVertex->set(xp, yp, zp, fTimeOffset );

  LogDebug("VertexGenerator") << "BeamProfileVtxGenerator: Vertex created "
                              << "at (" << xp << ", " << yp << ", "
                              << zp << ", " << fTimeOffset << ")";
  return fVertex;
}
BeamProfileVtxGenerator& BeamProfileVtxGenerator::operator= ( const BeamProfileVtxGenerator rhs) [private]

Copy assignment operator

void BeamProfileVtxGenerator::phi ( double  m = 0) [inline]

set phi in radian

Definition at line 51 of file BeamProfileVtxGenerator.h.

References fPhi, and m.

Referenced by BeamProfileVtxGenerator().

{fPhi=m;}
void BeamProfileVtxGenerator::psi ( double  m = 999) [inline]

set psi in radian

Definition at line 53 of file BeamProfileVtxGenerator.h.

References fPsi, and m.

Referenced by BeamProfileVtxGenerator().

{fPsi=m;}
void BeamProfileVtxGenerator::setType ( bool  m = true)

set type

Definition at line 202 of file BeamProfileVtxGenerator.cc.

References fRandom, fType, BaseEvtVtxGenerator::getEngine(), and alignCSCRings::s.

Referenced by BeamProfileVtxGenerator().

                                            { 

  fType = s;
  delete fRandom;
  
  if (fType == true)
    fRandom = new CLHEP::RandGaussQ(getEngine());
  else
    fRandom = new CLHEP::RandFlat(getEngine());
}
void BeamProfileVtxGenerator::sigmaX ( double  s = 1.0)

set resolution in X in cm

Definition at line 173 of file BeamProfileVtxGenerator.cc.

References fSigmaX, and alignCSCRings::s.

Referenced by BeamProfileVtxGenerator().

                                             { 

  if (s>=0) {
    fSigmaX = s; 
  } else {
    edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
                                       << " Illegal resolution in X " << s
                                       << "- set to default value 0 cm";
    fSigmaX = 0;
  }
}
void BeamProfileVtxGenerator::sigmaY ( double  s = 1.0)

set resolution in Y in cm

Definition at line 185 of file BeamProfileVtxGenerator.cc.

References fSigmaY, and alignCSCRings::s.

Referenced by BeamProfileVtxGenerator().

                                             { 

  if (s>=0) {
    fSigmaY = s; 
  } else {
    edm::LogWarning("VertexGenerator") << "Warning BeamProfileVtxGenerator:"
                                       << " Illegal resolution in Y " << s
                                       << "- set to default value 0 cm";
    fSigmaY = 0;
  }
}

Member Data Documentation

std::vector<double> BeamProfileVtxGenerator::fdistn [private]

Definition at line 71 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and newVertex().

Definition at line 65 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and eta().

Definition at line 69 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and newVertex().

Definition at line 64 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), meanX(), and newVertex().

Definition at line 64 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), meanY(), and newVertex().

Definition at line 64 of file BeamProfileVtxGenerator.h.

Referenced by beamPos(), BeamProfileVtxGenerator(), and newVertex().

Definition at line 65 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), newVertex(), and phi().

Definition at line 67 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), newVertex(), and psi().

CLHEP::HepRandom* BeamProfileVtxGenerator::fRandom [private]

Definition at line 72 of file BeamProfileVtxGenerator.h.

Referenced by newVertex(), setType(), and ~BeamProfileVtxGenerator().

Definition at line 63 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), newVertex(), and sigmaX().

Definition at line 63 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), newVertex(), and sigmaY().

Definition at line 65 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), eta(), and newVertex().

Reimplemented from BaseEvtVtxGenerator.

Definition at line 73 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and newVertex().

Definition at line 69 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), newVertex(), and setType().

Definition at line 70 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and newVertex().

Definition at line 70 of file BeamProfileVtxGenerator.h.

Referenced by BeamProfileVtxGenerator(), and newVertex().