CMS 3D CMS Logo

Public Member Functions

RPCRecHitStandardAlgo Class Reference

#include <RPCRecHitStandardAlgo.h>

Inheritance diagram for RPCRecHitStandardAlgo:
RPCRecHitBaseAlgo

List of all members.

Public Member Functions

virtual bool compute (const RPCRoll &roll, const RPCCluster &cluster, LocalPoint &point, LocalError &error) const
 standard local recHit computation
virtual bool compute (const RPCRoll &roll, const RPCCluster &cluster, const float &angle, const GlobalPoint &globPos, LocalPoint &point, LocalError &error) const
 RPCRecHitStandardAlgo (const edm::ParameterSet &config)
 Constructor.
virtual void setES (const edm::EventSetup &setup)
 Pass the Event Setup to the algo at each event.
virtual ~RPCRecHitStandardAlgo ()
 Destructor.

Detailed Description

Concrete implementation of RPCRecHitBaseAlgo.

Date:
2006/04/18 16:28:31
Revision:
1.1
Author:
M. Maggi -- INFN Bari

Definition at line 16 of file RPCRecHitStandardAlgo.h.


Constructor & Destructor Documentation

RPCRecHitStandardAlgo::RPCRecHitStandardAlgo ( const edm::ParameterSet config)

Constructor.

Definition at line 18 of file RPCRecHitStandardAlgo.cc.

                                                                          :
  RPCRecHitBaseAlgo(config) 
{
}
RPCRecHitStandardAlgo::~RPCRecHitStandardAlgo ( ) [virtual]

Destructor.

Definition at line 25 of file RPCRecHitStandardAlgo.cc.

{
}

Member Function Documentation

bool RPCRecHitStandardAlgo::compute ( const RPCRoll roll,
const RPCCluster cl,
LocalPoint Point,
LocalError error 
) const [virtual]

standard local recHit computation

Implements RPCRecHitBaseAlgo.

Definition at line 37 of file RPCRecHitStandardAlgo.cc.

References RPCRoll::centreOfStrip(), RPCCluster::firstStrip(), RPCCluster::lastStrip(), and RPCRoll::localError().

Referenced by compute().

{
  // Get Average Strip position
  float fstrip = (roll.centreOfStrip(cluster.firstStrip())).x();
  float lstrip = (roll.centreOfStrip(cluster.lastStrip())).x();
  float centreOfCluster = (fstrip + lstrip)/2;

  LocalPoint loctemp2(centreOfCluster,0.,0.);
 
  Point = loctemp2;
  error = roll.localError((cluster.firstStrip()+cluster.lastStrip())/2.);
  return true;
}
bool RPCRecHitStandardAlgo::compute ( const RPCRoll roll,
const RPCCluster cl,
const float &  angle,
const GlobalPoint globPos,
LocalPoint Point,
LocalError error 
) const [virtual]

local recHit computation accounting for track direction and absolute position

Implements RPCRecHitBaseAlgo.

Definition at line 55 of file RPCRecHitStandardAlgo.cc.

References compute(), and PV3DBase< T, PVType, FrameType >::z().

{

  // Glob Pos and angle not used so far...
  if (globPos.z()<0){ } // Fake use to avoid warnings
  if (angle<0.){ }      // Fake use to avoid warnings
  this->compute(roll,cl,Point,error);
  return true;
}
void RPCRecHitStandardAlgo::setES ( const edm::EventSetup setup) [virtual]

Pass the Event Setup to the algo at each event.

Implements RPCRecHitBaseAlgo.

Definition at line 31 of file RPCRecHitStandardAlgo.cc.

                                                            {
}