CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

CSCDriftSim Class Reference

#include <CSCDriftSim.h>

List of all members.

Public Member Functions

 CSCDriftSim ()
CSCDetectorHit getWireHit (const Local3DPoint &ionClusterPosition, const CSCLayer *, int wire, const PSimHit &simHit)
void setMagneticField (const MagneticField *field)
void setRandomEngine (CLHEP::HepRandomEngine &engine)
 ~CSCDriftSim ()

Private Member Functions

double avalancheCharge ()
double avgDrift () const
double avgDriftTimeHighB ()
double avgDriftTimeLowB ()
double avgPathLengthHighB ()
double avgPathLengthLowB ()
double driftSigma () const
double driftTimeSigmaHighB ()
double driftTimeSigmaLowB ()
double gasGain (const CSCDetId &id) const
double pathSigmaHighB ()
double pathSigmaLowB ()

Private Attributes

float bz
std::vector< double > dNdEIntegral
const double ELECTRON_DIFFUSION_COEFF
const double STEP_SIZE
const MagneticFieldtheMagneticField
CLHEP::RandFlat * theRandFlat
CLHEP::RandGaussQ * theRandGaussQ
double ycell
double zcell

Detailed Description

Simulate drift of electrons through CSC gas.

Author:
Rick Wilkinson

Someday this class will be subclassed into fast and slow sims, Now it's just slow, according to the rest of the ORCA developers...

Last mod:
30-Jun-00 ptc Doxygenate.
Bug-fix in ctor: first bin of integral was 'nan'.
Bug-trap in avalancheCharge(): very rarely could attempt access outside std::vector.

01-08/00 vin use binary search in avalancheCharge()

Definition at line 30 of file CSCDriftSim.h.


Constructor & Destructor Documentation

CSCDriftSim::CSCDriftSim ( )

Definition at line 20 of file CSCDriftSim.cc.

References dNdEIntegral, funct::exp(), i, N_INTEGRAL_STEPS, funct::pow(), and STEP_SIZE.

: bz(0.), // should make these local variables
  ycell(0.),
  zcell(0.),
  dNdEIntegral(N_INTEGRAL_STEPS, 0.),
  STEP_SIZE(0.01),
  ELECTRON_DIFFUSION_COEFF(0.0161),
  theMagneticField(0),
  theRandGaussQ(0),
  theRandFlat(0)
{
  // just initialize avalanche sim.  There has to be a better
  // way to take the integral of a function!
  double sum = 0.;
  int i;
  for(i = 0; i < N_INTEGRAL_STEPS; ++i) {
    if(i > 1) {
      double xx = STEP_SIZE * (double(i) - 0.5 );
      double dNdE = pow( xx, 0.38) * exp(-1.38*xx);

      sum += dNdE;
    }
    // store this value in the map
    dNdEIntegral[i] = sum;
  }

  // now normalize the whole map
  for(i =  0; i < N_INTEGRAL_STEPS; ++i) {
    dNdEIntegral[i] /= sum;
  }
}
CSCDriftSim::~CSCDriftSim ( )

Definition at line 53 of file CSCDriftSim.cc.

References theRandFlat, and theRandGaussQ.

{
  delete theRandGaussQ;
  delete theRandFlat;
}

Member Function Documentation

double CSCDriftSim::avalancheCharge ( ) [private]

Definition at line 144 of file CSCDriftSim.cc.

References dNdEIntegral, i, LogTrace, AlCaHLTBitMon_ParallelJobs::p, STEP_SIZE, theRandFlat, and x.

Referenced by getWireHit().

                                    {
  double returnVal = 0.;
  // pick a random value along the dNdE integral
  double x = theRandFlat->fire();
  size_t i;
  size_t isiz = dNdEIntegral.size();
  /*
  for(i = 0; i < isiz-1; ++i) {
    if(dNdEIntegral[i] > x) break;
  }
  */
  // return position of first element with a value >= x
  std::vector<double>::const_iterator p=lower_bound(dNdEIntegral.begin(),dNdEIntegral.end(),x);
  if (p==dNdEIntegral.end()) i=isiz-1;
  else i = p-dNdEIntegral.begin();
                              

  // now extrapolate between values
  if( i ==  isiz-1 ) {
    //edm::LogInfo("CSCDriftSim") << "Funky integral in CSCDriftSim " << x;
    returnVal = STEP_SIZE * double(i) * dNdEIntegral[i];
  }
  else {
    double x1 = dNdEIntegral[i];
    double x2 = dNdEIntegral[i+1];
    returnVal = STEP_SIZE * (double(i) + (x-x1)/(x2-x1)); 
  }
  LogTrace("CSCDriftSim") << "CSCDriftSim: avalanche fluc " << returnVal << "  " << x ;

  return returnVal;  
}
double CSCDriftSim::avgDrift ( ) const [private]

Definition at line 8 of file CSCDriftParam.cc.

References bz, and zcell.

Referenced by getWireHit().

                                  {
    /* Initialized data */

    static const double coeff[6] = { -.08618265,-.0020471805,-7.7771331e-4,
            -.0058433565,-5.0241491e-4,.0010223353 };


/* ! Parameterization of Lorentz drift in CSCs */
/* ***********************************************************************
 * DOC MC_ALL_BDFT                                                       *
 *                                                                       *
 * DOC  Function    : Parameterization of the Lorentz drift of electrons *
 * DOC                in the muon endcap CSCs.                           *
 *                                                                       *
 * DOC  References  : None                                               *
 *                                                                       *
 * DOC  Arguments   : BZ  - z componant of the local magnetic            *
 *                         field in kgauss                               *
 * DOC                ZCELL - distance from the anode wire in the wire   *
 *                           plane coordinate normalized to one          *
 * DOC  Errors      : None                                               *
 *                                                                       *
 * DOC  Returns     : Drift distance along the anode wire                *
 *                                                                       *
 * DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 * ***********************************************************************
 */




    const double x10 = 1.;
    const double x11 = bz / 40.;
    const double x12 = x11 * 2. * x11 - x10;
    const double x13 = x11 * 2. * x12 - x11;
    const double x14 = x11 * 2. * x13 - x12;
    const double x15 = x11 * 2. * x14 - x13;
    const double x16 = x11 * 2. * x15 - x14;
    const double x17 = x11 * 2. * x16 - x15;
    const double x18 = x11 * 2. * x17 - x16;
    const double x19 = x11 * 2. * x18 - x17;
    const double x20 = 1.;
    const double x21 = zcell;
    const double x22 = x21 * 2. * x21 - x20;
    const double x23 = x21 * 2. * x22 - x21;
    const double x24 = x21 * 2. * x23 - x22;
    const double x25 = x21 * 2. * x24 - x23;
    const double x26 = x21 * 2. * x25 - x24;
    const double x27 = x21 * 2. * x26 - x25;

  // not sure why we need a minus sign for avgDrift, but we think the
  // parametrization has a wrong sign.  For the +z endcap, bz should
  // be negative in local coord, so it should drift to the right
  // if drifting up (neg zcell), and to the left if it drifts down.  It doesn't.
  // the first-order term here is -0.08 * bz/40 * zcell

    return -1 * (coeff[0] * x11 * x21 + coeff[1] * x17 * x21 + coeff[2] * x21 +
            coeff[3] * x11 * x27 + coeff[4] * x22 + coeff[5] * x19);
}
double CSCDriftSim::avgDriftTimeHighB ( ) [private]

Definition at line 157 of file CSCDriftParamHighB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                      {
    /* Initialized data */

    static const double coeff[27] = { 22.384492,10.562894,14.032961,7.06233,
            3.5523289,-5.0176704,1.999075,1.0635552,-3.2770096,-2.7384958,
            .98411495,-2.0963696,-1.4006525,-.47542728,.64179451,-.80308436,
            .42964647,-.4153324,.50423068,.35049792,-.42595896,-.30947641,
            .16671267,-.21336584,.22979164,.23481052,.32550435 };

    /* System generated locals */
    float ret_val;

    /* Local variables */
    double x10, x11, x12, x13, x14, x15, x16, x17, x20, x21, x22,
            x23, x24, x25, x26, x27;

/* ! Parameterization of drift time - high field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_TIME                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift time                 *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/*                           (ycell=1 > d_acat)                         *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate (zcell=1 > d_anod/2.)     *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift time for high field CSC chambers             *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */


    x10 = 1.;
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * 2. * x11 - x10;
    x13 = x11 * 2. * x12 - x11;
    x14 = x11 * 2. * x13 - x12;
    x15 = x11 * 2. * x14 - x13;
    x16 = x11 * 2. * x15 - x14;
    x17 = x11 * 2. * x16 - x15;
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    x23 = x21 * 2. * x22 - x21;
    x24 = x21 * 2. * x23 - x22;
    x25 = x21 * 2. * x24 - x23;
    x26 = x21 * 2. * x25 - x24;
    x27 = x21 * 2. * x26 - x25;

    ret_val = coeff[0] + coeff[1] * x11 + coeff[2] * x21 + coeff[3] * x22 +
            coeff[4] * x23 + coeff[5] * x11 * x21 + coeff[6] * x24 + coeff[7]
            * x12 + coeff[8] * x11 * x22 + coeff[9] * x11 * x23 + coeff[10] *
            x25 + coeff[11] * x11 * x24 + coeff[12] * x11 * x25 + coeff[13] *
            x13 + coeff[14] * x12 * x21 + coeff[15] * x11 * x26 + coeff[16] *
            x26 + coeff[17] * x11 * x27 + coeff[18] * x17 * x21 + coeff[19] *
            x15 * x21 + coeff[20] * x12 * x22 + coeff[21] * x12 * x23 + coeff[
            22] * x27 + coeff[23] * x14 * x22 + coeff[24] * x16 * x21 + coeff[
            25] * x17 + coeff[26] * x17 * x22;

    return ret_val;
}
double CSCDriftSim::avgDriftTimeLowB ( ) [private]

Definition at line 155 of file CSCDriftParamLowB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                     {
    /* Initialized data */

    static const double coeff[20] = { 42.981588,25.732805,26.539129,16.719016,
            10.862044,7.4859085,5.0353142,3.3620548,1.9057762,2.2207695,
            -2.6787582,1.2977292,-.8358091,1.2452612,.74813469,-.57581194,
            .32705275,-.85426489,-.55688158,-.38384903 };

    /* System generated locals */
    float ret_val;

    /* Local variables */
    double x10, x11, x12, x13, x14, x20, x21, x22, x23, x24, x25,
            x26, x27, x28, x29;

/* ! Parameterization of drift time - low field chambers */
/* ***********************************************************************
 */
/* DOC MC_BLOW_TIME                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift time                 *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/*                           (ycell=1 > d_acat)                         *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate (zcell=1 > d_anod/2.)     *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift time for low field CSC chambers              *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */




    x10 = 1.;
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * 2. * x11 - x10;
    x13 = x11 * 2. * x12 - x11;
    x14 = x11 * 2. * x13 - x12;
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    x23 = x21 * 2. * x22 - x21;
    x24 = x21 * 2. * x23 - x22;
    x25 = x21 * 2. * x24 - x23;
    x26 = x21 * 2. * x25 - x24;
    x27 = x21 * 2. * x26 - x25;
    x28 = x21 * 2. * x27 - x26;
    x29 = x21 * 2. * x28 - x27;

    ret_val = coeff[0] + coeff[1] * x11 + coeff[2] * x21 + coeff[3] * x22 +
            coeff[4] * x23 + coeff[5] * x24 + coeff[6] * x25 + coeff[7] * x26
            + coeff[8] * x12 + coeff[9] * x27 + coeff[10] * x11 * x21 + coeff[
            11] * x28 + coeff[12] * x13 + coeff[13] * x12 * x21 + coeff[14] *
            x29 + coeff[15] * x13 * x21 + coeff[16] * x14 + coeff[17] * x11 *
            x22 + coeff[18] * x11 * x23 + coeff[19] * x11 * x24;
    return ret_val;
}
double CSCDriftSim::avgPathLengthHighB ( ) [private]

Definition at line 4 of file CSCDriftParamHighB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                       {
    /* Initialized data */

    static const double coeff[18] = { .16916627,.11057547,.054287448,.01179527,
            .0062073273,-.013570915,-.0027121772,-.0053792764,-.0027452986,
            -.0020556715,.0021511659,.0011376412,.0026183373,.0017980602,
            -.0012975418,-.0010798782,-.0012322628,-8.3635924e-4 };

    /* Local variables */
    double x10, x11, x12, x13, x14, x20, x21, x22, x23, x24, x25;

/* ! Parameterization of drift path length - high field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_SLEN                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift path length          *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate                           *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift path length for high field CSC chambers      *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */



    x10 = 1.;
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * 2. * x11 - x10;
    x13 = x11 * 2. * x12 - x11;
    x14 = x11 * 2. * x13 - x12;
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    x23 = x21 * 2. * x22 - x21;
    x24 = x21 * 2. * x23 - x22;
    x25 = x21 * 2. * x24 - x23;

    return coeff[0] + coeff[1] * x11 + coeff[2] * x21 + coeff[3] * x22 +
            coeff[4] * x12 + coeff[5] * x11 * x21 + coeff[6] * x13 + coeff[7]
            * x12 * x22 + coeff[8] * x12 * x23 + coeff[9] * x11 * x24 + coeff[
            10] * x12 * x21 + coeff[11] * x14 + coeff[12] * x11 * x22 + coeff[
            13] * x13 * x22 + coeff[14] * x13 * x21 + coeff[15] * x12 * x24 +
            coeff[16] * x11 * x25 + coeff[17] * x11 * x23;
}
double CSCDriftSim::avgPathLengthLowB ( ) [private]

Definition at line 7 of file CSCDriftParamLowB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                      {

    static  const double coeff[7] = { .28872209,.21705601,.063908389,-.012224924,
            .012901814,-.0089355058,-.015769145 };

    // System generated locals
    float ret_val;

    // Local variables
    double x10, x11, x12, x13, x20, x21, x22;

/* ! Parameterization of drift path length - low field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_SLEN                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift path length          *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate                           *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift path length for low field CSC chambers      *
*/
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */


    x10 = 1.;
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * 2. * x11 - x10;
    x13 = x11 * 2. * x12 - x11;
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    ret_val = coeff[0] + coeff[1] * x11 + coeff[2] * x21 + coeff[3] * x12 *
            x22 + coeff[4] * x22 + coeff[5] * x13 + coeff[6] * x11 * x21;
    return ret_val;
}
double CSCDriftSim::driftSigma ( ) const [private]

Definition at line 69 of file CSCDriftParam.cc.

References bz, and zcell.

Referenced by getWireHit().

                                     {
    // Initialized data 

    static const double coeff[3] = { .01069525,.010364504,.0021662697 };


/* ! Parameterization of Lorentz drift error in CSCs */
/*************************************************************************
 * DOC MC_ALL_BSIG                                                       *
 *                                                                       *
 * DOC  Function    : Parameterization of the dispersion in Lorentz Drift*
 * DOC                in the muon endcap CSCs.                           *
 *                                                                       *
 * DOC  References  : None                                               *
 *                                                                       *
 * DOC  Arguments   : BZ  - z componant of the local magnetic            *
 * DOC                ZCELL - distance from the anode wire in the wire   *
 *                           plane coordinate normalized to 1            *
 * DOC  Errors      : None                                               *
 *                                                                       *
 * DOC  Returns     : Dispersion in Lorentz drift along anode wire       *
 *                                                                       *
 * DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 * ***********************************************************************
 */


    const double x11 = bz / 40.;
    const double x20 = 1.;
    const double x21 = zcell;
    const double x22 = x21 * 2. * x21 - x20;

    return coeff[0] + coeff[1] * x11 * x21 + coeff[2] * x22;
} 
double CSCDriftSim::driftTimeSigmaHighB ( ) [private]

Definition at line 245 of file CSCDriftParamHighB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                        {
    /* Initialized data */

    static const double coeff[17] = { 5.5533465,3.3733352,3.776603,2.2673355,
            1.3401485,.84209333,-.71621378,.57572407,-.52313936,-.78790514,
            -.71786066,.43370011,.29223306,-.37791975,.21121024,.31513644,
            .25382701 };

    /* System generated locals */
    float ret_val;

    /* Local variables */
    double x10, x11, x12, x13, x14, x15, x16, x17, x18, /*x19,*/ x20,
            x21, x22, x23, x24, x25, x26, x27, x28, x29;

/* ! Parameterization of drift time dispersion- high field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_TSIG                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift time dispersion      *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/*                           (ycell=1 > d_acat)                         *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate (zcell=1 > d_anod/2.)     *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift time dispersion for high field CSC chambers  *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */



    x10 = 1.;
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * 2. * x11 - x10;
    x13 = x11 * 2. * x12 - x11;
    x14 = x11 * 2. * x13 - x12;
    x15 = x11 * 2. * x14 - x13;
    x16 = x11 * 2. * x15 - x14;
    x17 = x11 * 2. * x16 - x15;
    x18 = x11 * 2. * x17 - x16;
    //x19 = x11 * 2. * x18 - x17; //not used later
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    x23 = x21 * 2. * x22 - x21;
    x24 = x21 * 2. * x23 - x22;
    x25 = x21 * 2. * x24 - x23;
    x26 = x21 * 2. * x25 - x24;
    x27 = x21 * 2. * x26 - x25;
    x28 = x21 * 2. * x27 - x26;
    x29 = x21 * 2. * x28 - x27;

    ret_val = coeff[0] * x21 + coeff[1] + coeff[2] * x22 + coeff[3] * x23 +
            coeff[4] * x24 + coeff[5] * x25 + coeff[6] * x11 * x23 + coeff[7]
            * x26 + coeff[8] * x11 * x25 + coeff[9] * x11 * x24 + coeff[10] *
            x11 * x22 + coeff[11] * x27 + coeff[12] * x28 + coeff[13] * x11 *
            x26 + coeff[14] * x29 + coeff[15] * x16 * x21 + coeff[16] * x18 *
            x21;

    return ret_val;
} 
double CSCDriftSim::driftTimeSigmaLowB ( ) [private]

Definition at line 239 of file CSCDriftParamLowB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                       {
    /* Initialized data */

    static  const double coeff[19] = { 6.2681223,3.5916437,5.5425809,4.6974052,
            3.8472392,3.1019155,2.4323913,1.8695623,1.3909068,.99056625,
            .52066638,.64671229,.68023389,.40512251,.40202738,.23908531,
            -.41245784,-.32196924,-.29890696 };

    /* Local variables */
    double /*x10,*/ x11, x20, x21, x22, x23, x24, x25, x26, x27, x28,
      x29, x210, x211, x212;

/* ! Parameterization of drift time dispersion- low field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_TSIG                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift time dispersion      *
 */
/* DOC                in the muon endcap CSCs.                           *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/*                           (ycell=1 > d_acat)                         *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate (zcell=1 > d_anod/2.)     *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Drift time dispersion for low field CSC chambers   *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */



    //x10 = 1.; //not used
    x11 = fabs(ycell) * 2. - 1.;
    x20 = 1.;
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * 2. * x21 - x20;
    x23 = x21 * 2. * x22 - x21;
    x24 = x21 * 2. * x23 - x22;
    x25 = x21 * 2. * x24 - x23;
    x26 = x21 * 2. * x25 - x24;
    x27 = x21 * 2. * x26 - x25;
    x28 = x21 * 2. * x27 - x26;
    x29 = x21 * 2. * x28 - x27;
    x210 = x21 * 2. * x29 - x28;
    x211 = x21 * 2. * x210 - x29;
    x212 = x21 * 2. * x211 - x210;

    return coeff[0] * x21 + coeff[1] + coeff[2] * x22 + coeff[3] * x23 +
            coeff[4] * x24 + coeff[5] * x25 + coeff[6] * x26 + coeff[7] * x27
            + coeff[8] * x28 + coeff[9] * x29 + coeff[10] * x11 + coeff[11] *
            x11 * x21 + coeff[12] * x210 + coeff[13] * x11 * x22 + coeff[14] *
             x211 + coeff[15] * x212 + coeff[16] * x11 * x25 + coeff[17] *
            x11 * x27 + coeff[18] * x11 * x26;
}
double CSCDriftSim::gasGain ( const CSCDetId id) const [private]

Definition at line 177 of file CSCDriftSim.cc.

References query::result, CSCDetId::ring(), relativeConstraints::ring, and CSCDetId::station().

Referenced by getWireHit().

{
  double result = 130000.;  // 1.30 E05
  // if ME1/1, add some extra gas gain to compensate
  // for a smaller gas gap
  int ring = detId.ring();
  if(detId.station() == 1 && (ring == 1 || ring == 4))
  {
    result = 213000.; // 2.13 E05 ( to match real world as of Jan-2011)
  }
  return result;
}
CSCDetectorHit CSCDriftSim::getWireHit ( const Local3DPoint ionClusterPosition,
const CSCLayer layer,
int  wire,
const PSimHit simHit 
)

takes a point, and creates a signal on the wire

Definition at line 68 of file CSCDriftSim.cc.

References CSCChamberSpecs::anodeCathodeSpacing(), avalancheCharge(), avgDrift(), avgDriftTimeHighB(), avgDriftTimeLowB(), avgPathLengthHighB(), avgPathLengthLowB(), bz, CSCLayer::chamber(), DeDxDiscriminatorTools::charge(), driftSigma(), driftTimeSigmaHighB(), driftTimeSigmaLowB(), e_SI, ELECTRON_DIFFUSION_COEFF, gasGain(), relativeConstraints::geom, CSCLayer::geometry(), CSCLayer::id(), MagneticField::inKGauss(), LogTrace, max(), pathSigmaHighB(), pathSigmaLowB(), pos, idealTransformation::rotation, CSCChamber::specs(), mathSSE::sqrt(), GeomDet::surface(), lumiQTWidget::t, theMagneticField, theRandGaussQ, PSimHit::tof(), Surface::toGlobal(), GeomDet::toLocal(), CSCLayerGeometry::wireAngle(), CSCChamberSpecs::wireSpacing(), PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), ycell, CSCLayerGeometry::yOfWire(), PV3DBase< T, PVType, FrameType >::z(), z, and zcell.

Referenced by CSCWireHitSim::simulate().

                                                                   {

  const CSCChamberSpecs * specs = layer->chamber()->specs();
  const CSCLayerGeometry * geom = layer->geometry();
  math::LocalPoint clusterPos(pos.x(), pos.y(), pos.z());
  LogTrace("CSCDriftSim") << "CSCDriftSim: ionization cluster at: " <<  pos; 
  // set the coordinate system with the x-axis along the nearest wire,
  // with the origin in the center of the chamber, on that wire.
  math::LocalVector yShift(0, -1.*geom->yOfWire(nearestWire), 0.);
  ROOT::Math::RotationZ rotation(-1.*geom->wireAngle());

  clusterPos = yShift + clusterPos;
  clusterPos = rotation * clusterPos;
  GlobalPoint globalPosition = layer->surface().toGlobal(pos);
  assert(theMagneticField != 0);
 
  //  bz = theMagneticField->inTesla(globalPosition).z() * 10.;

  // We need magnetic field in _local_ coordinates
  // Interface now allows access in kGauss directly.
  bz = layer->toLocal( theMagneticField->inKGauss( globalPosition ) ).z();

  // these subroutines label the coordinates in GEANT coords...
  ycell = clusterPos.z() / specs->anodeCathodeSpacing();
  zcell = 2.*clusterPos.y() / specs->wireSpacing();

  LogTrace("CSCDriftSim") << "CSCDriftSim: bz " << bz <<" avgDrift " << avgDrift()
       << " wireAngle " << geom->wireAngle()
       << " ycell " << ycell << " zcell " << zcell;

  double avgPathLength, pathSigma, avgDriftTime, driftTimeSigma;
  static const float B_FIELD_CUT = 15.f;
  if(fabs(bz) < B_FIELD_CUT) {
    avgPathLength  = avgPathLengthLowB();
    pathSigma      = pathSigmaLowB();
    avgDriftTime   = avgDriftTimeLowB();
    driftTimeSigma = driftTimeSigmaLowB();
  }
  else {
    avgPathLength  = avgPathLengthHighB();
    pathSigma      = pathSigmaHighB();
    avgDriftTime   = avgDriftTimeHighB();
    driftTimeSigma = driftTimeSigmaHighB();
  }

  // electron drift path length 
  double pathLength = std::max(theRandGaussQ->fire(avgPathLength, pathSigma), 0.);

  // electron drift distance along the anode wire, including diffusion
  double diffusionSigma = ELECTRON_DIFFUSION_COEFF * sqrt(pathLength);
  double x = clusterPos.x() + theRandGaussQ->fire(avgDrift(), driftSigma())
                            + theRandGaussQ->fire(0., diffusionSigma);

  // electron drift time
  double driftTime  = std::max(theRandGaussQ->fire(avgDriftTime, driftTimeSigma), 0.);

  //@@ Parameters which should be defined outside the code
  // f_att is the fraction of drift electrons lost due to attachment
  //static const double f_att = 0.5;
  static const double f_collected = 0.82;

  // Avalanche charge, with fluctuation ('avalancheCharge()' is the fluctuation generator!)
  //double charge = avalancheCharge() * f_att * f_collected * gasGain(layer->id()) * e_SI * 1.e15;
  // doing fattachment by random chance of killing
  double charge = avalancheCharge() * f_collected * gasGain(layer->id()) * e_SI * 1.e15;

  float t = simHit.tof() + driftTime;
  LogTrace("CSCDriftSim") << "CSCDriftSim: tof = " << simHit.tof() << 
      " driftTime = " << driftTime <<
      " MEDH = " << CSCDetectorHit(nearestWire, charge, x, t, &simHit);
  return CSCDetectorHit(nearestWire, charge, x, t, &simHit);
}
double CSCDriftSim::pathSigmaHighB ( ) [private]

Definition at line 76 of file CSCDriftParamHighB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                   {
    /* Initialized data */

    static const double coeff[9] = { .0049089564,.0091482062,.0024036507,
            .0065285652,.0041487742,-.0038102526,-.0043923587,.0019230151,
            .0013543258 };

    /* System generated locals */
    float ret_val;

    /* Local variables */
    double /*x10,*/ x11, x12, x13, x14, x15, x16, /*x20,*/ x21, x22, x23,
            x24, x25, x26, x27, x28, x29;

/* ! Parameterization of path length dispersion- high field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_SSIG                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift path length          *
 */
/* DOC                dispersion in the muon endcap CSCs.                *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate                           *
*/
/*           **NOTE** Both distances normalize to cell dim=1x1          *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Path length dispersion for high field CSC chambers *
 */
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */


    //x10 = 1.; //not used later
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * x11;
    x13 = x11 * x12;
    x14 = x11 * x13;
    x15 = x11 * x14;
    x16 = x11 * x15;
    //x20 = 1.; //not used later
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * x21;
    x23 = x21 * x22;
    x24 = x21 * x23;
    x25 = x21 * x24;
    x26 = x21 * x25;
    x27 = x21 * x26;
    x28 = x21 * x27;
    x29 = x21 * x28;

    ret_val = coeff[0] + coeff[1] * x21 + coeff[2] * x11 + coeff[3] * x22 +
            coeff[4] * x11 * x21 + coeff[5] * x16 * x22 + coeff[6] * x16 *
            x23 + coeff[7] * x11 * x22 + coeff[8] * x29;

    return ret_val;
}
double CSCDriftSim::pathSigmaLowB ( ) [private]

Definition at line 70 of file CSCDriftParamLowB.cc.

References ycell, and zcell.

Referenced by getWireHit().

                                  {
    /* Initialized data */

    static  const double coeff[12] = { .002854441,8.701339e-4,.0053064193,
            .0012356508,9.8627318e-4,.0013802449,-5.4633755e-4,.0026078648,
            .0044171026,6.2029063e-4,.0011392474,-.0056275595 };

    /* System generated locals */
    float ret_val;

    /* Local variables */
    double /*x10,*/ x11, x12, x13, x14, x15, x16, /*x20,*/ x21, x22, x23,
            x24, x25, x26, x27, x28, x29, x210, x211, x212, x213;

/* ! Parameterization of path length dispersion- low field chambers */
/* ***********************************************************************
 */
/* DOC MC_BHGH_SSIG                                                      *
 */
/*                                                                      *
*/
/* DOC  Function    : Parameterization of the drift path length          *
 */
/* DOC                dispersion in the muon endcap CSCs.                *
 */
/*                                                                      *
*/
/* DOC  References  : None                                               *
 */
/*                                                                      *
*/
/* DOC  Arguments   : YCELL - distance from the anode wire in the        *
 */
/*                           anode-cathode coordinate plane             *
*/
/* DOC                ZCELL - distance from the anode wire in the wire   *
 */
/*                           plane coordinate                           *
*/
/*           **NOTE** Both distances normalize to cell dim=1x1          *
*/
/* DOC  Errors      : None                                               *
 */
/*                                                                      *
*/
/* DOC  Returns     : Path length dispersion for low field CSC chambers *
*/
/*                                                                      *
*/
/* DOC  Created     : 15-OCT-1996   Author : Jeff Rowe                   *
 */
/* ***********************************************************************
 */


    //x10 = 1.; //not used later
    x11 = fabs(ycell) * 2. - 1.;
    x12 = x11 * x11;
    x13 = x11 * x12;
    x14 = x11 * x13;
    x15 = x11 * x14;
    x16 = x11 * x15;
    //x20 = 1.; //not used later
    x21 = fabs(zcell) * 2. - 1.;
    x22 = x21 * x21;
    x23 = x21 * x22;
    x24 = x21 * x23;
    x25 = x21 * x24;
    x26 = x21 * x25;
    x27 = x21 * x26;
    x28 = x21 * x27;
    x29 = x21 * x28;
    x210 = x21 * x29;
    x211 = x21 * x210;
    x212 = x21 * x211;
    x213 = x21 * x212;
    ret_val = coeff[0] + coeff[1] * x23 + coeff[2] * x24 + coeff[3] * x13 +
            coeff[4] * x11 * x21 + coeff[5] * x11 * x22 + coeff[6] * x16 +
            coeff[7] * x213 + coeff[8] * x212 + coeff[9] * x21 + coeff[10] *
            x11 * x23 + coeff[11] * x26;

    return ret_val;
}
void CSCDriftSim::setMagneticField ( const MagneticField field) [inline]

Definition at line 43 of file CSCDriftSim.h.

References theMagneticField.

Referenced by CSCDigitizer::setMagneticField().

{theMagneticField = field;}
void CSCDriftSim::setRandomEngine ( CLHEP::HepRandomEngine &  engine)

Definition at line 60 of file CSCDriftSim.cc.

References theRandFlat, and theRandGaussQ.

Referenced by CSCWireHitSim::setRandomEngine().

{
  theRandGaussQ = new CLHEP::RandGaussQ(engine);
  theRandFlat = new CLHEP::RandFlat(engine);
}

Member Data Documentation

float CSCDriftSim::bz [private]

Definition at line 63 of file CSCDriftSim.h.

Referenced by avgDrift(), driftSigma(), and getWireHit().

std::vector<double> CSCDriftSim::dNdEIntegral [private]

Definition at line 67 of file CSCDriftSim.h.

Referenced by avalancheCharge(), and CSCDriftSim().

const double CSCDriftSim::ELECTRON_DIFFUSION_COEFF [private]

Definition at line 70 of file CSCDriftSim.h.

Referenced by getWireHit().

const double CSCDriftSim::STEP_SIZE [private]

Definition at line 68 of file CSCDriftSim.h.

Referenced by avalancheCharge(), and CSCDriftSim().

Definition at line 72 of file CSCDriftSim.h.

Referenced by getWireHit(), and setMagneticField().

CLHEP::RandFlat* CSCDriftSim::theRandFlat [private]

Definition at line 75 of file CSCDriftSim.h.

Referenced by avalancheCharge(), setRandomEngine(), and ~CSCDriftSim().

CLHEP::RandGaussQ* CSCDriftSim::theRandGaussQ [private]

Definition at line 74 of file CSCDriftSim.h.

Referenced by getWireHit(), setRandomEngine(), and ~CSCDriftSim().

double CSCDriftSim::ycell [private]
double CSCDriftSim::zcell [private]