CMS 3D CMS Logo

Public Member Functions | Private Attributes

MuRingForwardLayer Class Reference

#include <MuRingForwardLayer.h>

Inheritance diagram for MuRingForwardLayer:
RingedForwardLayer ForwardDetLayer DetLayer GeometricSearchDet

List of all members.

Public Member Functions

virtual const std::vector
< const GeomDet * > & 
basicComponents () const
virtual std::vector< DetWithStatecompatibleDets (const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
virtual const std::vector
< const GeometricSearchDet * > & 
components () const
 Returns basic components, if any.
virtual std::vector< DetGroupgroupedCompatibleDets (const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
virtual bool hasGroups () const
 MuRingForwardLayer (const std::vector< const ForwardDetRing * > &rings)
 Constructor, takes ownership of pointers.
virtual const std::vector
< const ForwardDetRing * > & 
rings () const
 Return the vector of rings.
virtual SubDetector subDetector () const
 The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
virtual ~MuRingForwardLayer ()

Private Attributes

bool isOverlapping
std::vector< const GeomDet * > theBasicComps
BaseBinFinder< double > * theBinFinder
std::vector< const
GeometricSearchDet * > 
theComponents
std::vector< const
ForwardDetRing * > 
theRings

Detailed Description

A plane composed of disks (MuRingForwardDisk). Represents forward muon CSC/RPC stations.

Date:
2007/06/14 17:22:41
Revision:
1.9
Author:
N. Amapane - INFN Torino

Definition at line 20 of file MuRingForwardLayer.h.


Constructor & Destructor Documentation

MuRingForwardLayer::MuRingForwardLayer ( const std::vector< const ForwardDetRing * > &  rings)

Constructor, takes ownership of pointers.

Definition at line 26 of file MuRingForwardLayer.cc.

References basicComponents(), isOverlapping, RBorderFinder::isROverlapping(), RBorderFinder::isRPeriodic(), LogTrace, max(), metname, min, pos, makeMuonMisalignmentScenario::rot, ForwardDetLayer::setSurface(), ForwardDetLayer::specificSurface(), AlCaHLTBitMon_QueryRunRegistry::string, theBasicComps, theBinFinder, theRings, and zPos.

                                                                                 :
  theRings(rings),
  theComponents(theRings.begin(),theRings.end()),
  theBinFinder(0),
  isOverlapping(false) 
{

  const std::string metname = "Muon|RecoMuon|RecoMuonDetLayers|MuRingForwardLayer";

  // Initial values for R and Z bounds
  float theRmin = rings.front()->basicComponents().front()->position().perp(); 
  float theRmax = theRmin;
  float theZmin = rings.front()->position().z();
  float theZmax = theZmin;

  // Cache chamber pointers (the basic components_)
  // and find extension in R and Z
  for (vector<const ForwardDetRing*>::const_iterator it=rings.begin();
       it!=rings.end(); it++) {
    vector<const GeomDet*> tmp2 = (*it)->basicComponents();
    theBasicComps.insert(theBasicComps.end(),tmp2.begin(),tmp2.end());
    
    theRmin = min( theRmin, (*it)->specificSurface().innerRadius());
    theRmax = max( theRmax, (*it)->specificSurface().outerRadius());
    float halfThick = (*it)->surface().bounds().thickness()/2.;
    float zCenter = (*it)->surface().position().z();
    theZmin = min( theZmin, zCenter-halfThick);
    theZmax = max( theZmax, zCenter+halfThick); 
  }  
  
  RBorderFinder bf(theRings);
  isOverlapping = bf.isROverlapping();
  theBinFinder = new GeneralBinFinderInR<double>(bf);

  // Build surface
  
  float zPos = (theZmax+theZmin)/2.;
  PositionType pos(0.,0.,zPos);
  RotationType rot;

  setSurface(new BoundDisk( pos, rot, 
                            SimpleDiskBounds( theRmin, theRmax, 
                                              theZmin-zPos, theZmax-zPos)));


   
  LogTrace(metname) << "Constructing MuRingForwardLayer: "
                    << basicComponents().size() << " Dets " 
                    << theRings.size() << " Rings "
                    << " Z: " << specificSurface().position().z()
                    << " R1: " << specificSurface().innerRadius()
                    << " R2: " << specificSurface().outerRadius()
                    << " Per.: " << bf.isRPeriodic()
                    << " Overl.: " << bf.isROverlapping();
}
MuRingForwardLayer::~MuRingForwardLayer ( ) [virtual]

Definition at line 83 of file MuRingForwardLayer.cc.

References i, theBinFinder, and theRings.

                                       {
  delete theBinFinder;
  for (vector <const ForwardDetRing*>::iterator i = theRings.begin();
       i<theRings.end(); i++) {delete *i;}
}

Member Function Documentation

virtual const std::vector<const GeomDet*>& MuRingForwardLayer::basicComponents ( ) const [inline, virtual]

Implements GeometricSearchDet.

Definition at line 32 of file MuRingForwardLayer.h.

References theBasicComps.

Referenced by MuRingForwardLayer(), and MuRingForwardDoubleLayer::selfTest().

{return theBasicComps;}
vector< GeometricSearchDet::DetWithState > MuRingForwardLayer::compatibleDets ( const TrajectoryStateOnSurface startingState,
const Propagator prop,
const MeasurementEstimator est 
) const [virtual]

Returns all Dets compatible with a trajectory state according to the estimator est. The startingState should be propagated to the surface of each compatible Det using the Propagator passed as an argument. The default implementation should be overridden in dets with specific surface types to avoid propagation to a generic Surface

Reimplemented from GeometricSearchDet.

Definition at line 91 of file MuRingForwardLayer.cc.

References BaseBinFinder< T >::binIndex(), ForwardDetLayer::compatible(), GeometricSearchDet::compatibleDets(), TrajectoryStateOnSurface::globalPosition(), TrajectoryStateOnSurface::hasError(), TrajectoryStateOnSurface::localError(), LogTrace, metname, PV3DBase< T, PVType, FrameType >::perp(), LocalTrajectoryError::positionError(), query::result, ForwardDetLayer::specificSurface(), ForwardDetRing::specificSurface(), mathSSE::sqrt(), AlCaHLTBitMon_QueryRunRegistry::string, ForwardDetLayer::surface(), theBinFinder, theRings, toLocal(), LocalError::xx(), and LocalError::yy().

Referenced by MuRingForwardDoubleLayer::groupedCompatibleDets().

                                                                          {
  
  const std::string metname = "Muon|RecoMuon|RecoMuonDetLayers|MuRingForwardLayer";
  vector<DetWithState> result; 
  
  
  LogTrace(metname) << "MuRingForwardLayer::compatibleDets," 
                    << " R1 " << specificSurface().innerRadius()
                    << " R2: " << specificSurface().outerRadius()
                    << " FTS at R: " << startingState.globalPosition().perp();
  
  pair<bool, TrajectoryStateOnSurface> compat =
    compatible(startingState, prop, est);
  
  if (!compat.first) {
    
    LogTrace(metname) << "     MuRingForwardLayer::compatibleDets: not compatible"
                      << " (should not have been selected!)";
    return result;
  }
  
  
  TrajectoryStateOnSurface& tsos = compat.second;
  
  int closest = theBinFinder->binIndex(tsos.globalPosition().perp());
  const ForwardDetRing* closestRing = theRings[closest];
  
  // Check the closest ring
  
  LogTrace(metname) << "     MuRingForwardLayer::fastCompatibleDets, closestRing: "
                    << closest
                    << " R1 " << closestRing->specificSurface().innerRadius()
                    << " R2: " << closestRing->specificSurface().outerRadius()
                    << " FTS R: " << tsos.globalPosition().perp();
  if (tsos.hasError()) {
    LogTrace(metname)  << " sR: " << sqrt(tsos.localError().positionError().yy())
                       << " sX: " << sqrt(tsos.localError().positionError().xx());
  }
  LogTrace(metname) << endl;
   
   
  result = closestRing->compatibleDets(tsos, prop, est);
   
  int nclosest = result.size(); int nnextdet=0; // MDEBUG counters
   
  //FIXME: if closest is not compatible next cannot be either?
   
  // Use state on layer surface. Note that local coordinates and errors
  // are the same on the layer and on all rings surfaces, since 
  // all BoundDisks are centered in 0,0 and have the same rotation.
  // CAVEAT: if the rings are not at the same Z, the local position and error
  // will be "Z-projected" to the rings. This is a fairly good approximation.
  // However in this case additional propagation will be done when calling
  // compatibleDets.
  GlobalPoint startPos = tsos.globalPosition();
  LocalPoint nextPos(surface().toLocal(startPos));
   
  for (unsigned int idet=closest+1; idet < theRings.size(); idet++) {
    bool inside = false;
    if (tsos.hasError()) {
      inside=theRings[idet]->specificSurface().bounds().inside(nextPos,tsos.localError().positionError());
    } else {
      inside=theRings[idet]->specificSurface().bounds().inside(nextPos);
    }
    if (inside){
      LogTrace(metname) << "     MuRingForwardLayer::fastCompatibleDets:NextRing" << idet
                        << " R1 " << theRings[idet]->specificSurface().innerRadius()
                        << " R2: " << theRings[idet]->specificSurface().outerRadius()
                        << " FTS R " << nextPos.perp();
      nnextdet++;      
      vector<DetWithState> nextRodDets =
        theRings[idet]->compatibleDets(tsos, prop, est);
      if (nextRodDets.size()!=0) {
        result.insert( result.end(), 
                       nextRodDets.begin(), nextRodDets.end());
      } else {
        break;
      }
    }
  }
   
  for (int idet=closest-1; idet >= 0; idet--) {
    bool inside = false;
    if (tsos.hasError()) {
      inside=theRings[idet]->specificSurface().bounds().inside(nextPos,tsos.localError().positionError());
    } else {
      inside=theRings[idet]->specificSurface().bounds().inside(nextPos);
    }
    if (inside){
      LogTrace(metname) << "     MuRingForwardLayer::fastCompatibleDets:PreviousRing:" << idet
                        << " R1 " << theRings[idet]->specificSurface().innerRadius()
                        << " R2: " << theRings[idet]->specificSurface().outerRadius()
                        << " FTS R " << nextPos.perp();
      nnextdet++;
      vector<DetWithState> nextRodDets =
        theRings[idet]->compatibleDets(tsos, prop, est);
      if (nextRodDets.size()!=0) {
        result.insert( result.end(), 
                       nextRodDets.begin(), nextRodDets.end());
      } else {
        break;
      }
    }
  }
   
  LogTrace(metname) << "     MuRingForwardLayer::fastCompatibleDets: found: "
                    << result.size()
                    << " on closest: " << nclosest
                    << " # checked rings: " << 1 + nnextdet;
   
  return result;
}
const vector< const GeometricSearchDet * > & MuRingForwardLayer::components ( ) const [virtual]

Returns basic components, if any.

Returns direct components, if any

Implements GeometricSearchDet.

Definition at line 228 of file MuRingForwardLayer.cc.

References theComponents.

                                     {
  return theComponents;
}
vector< DetGroup > MuRingForwardLayer::groupedCompatibleDets ( const TrajectoryStateOnSurface startingState,
const Propagator prop,
const MeasurementEstimator est 
) const [virtual]

Similar to compatibleDets(), but the compatible Dets are grouped in one or more groups. Dets are put in the same group if they are mutually exclusive for track crossing, i.e. a reconstructible track cannot cross more than one Det from a group. Pathological tracks (spirals etc.) can of course violate this rule.
The DetGroups are sorted in the sequence of crossing by a track. In order to define the direction of crossing the Propagator used in this method should have a defined direction() : either "alongMomentum" or "oppositeToMomentum" but not "anyDirection".
The three signatures of this method differ by the input trajectory state arguments: the starting state can be a TrajectoryStateOnSurface or a FreeTrajectoryState, and the state on this CompositeDet may be already known or not. The last two arguments are as for the method compatibleDets().
First signature: The first argument is a TrajectoryStateOnSurface, usually not on the surface of this CompositeDet.

Reimplemented from GeometricSearchDet.

Definition at line 208 of file MuRingForwardLayer.cc.

References gather_cfg::cout.

                                                                                  {
  // FIXME should return only 1 group 
  cout << "dummy implementation of MuRingForwardLayer::groupedCompatibleDets()" << endl;
  return vector<DetGroup>();
}
bool MuRingForwardLayer::hasGroups ( ) const [virtual]

Implements GeometricSearchDet.

Definition at line 217 of file MuRingForwardLayer.cc.

                                         {
  // FIXME : depending on isOverlapping?
  return false;
}
virtual const std::vector<const ForwardDetRing*>& MuRingForwardLayer::rings ( ) const [inline, virtual]

Return the vector of rings.

Definition at line 57 of file MuRingForwardLayer.h.

References theRings.

Referenced by MuRingForwardDoubleLayer::isCrack().

{return theRings;}
GeomDetEnumerators::SubDetector MuRingForwardLayer::subDetector ( ) const [virtual]

The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)

Implements DetLayer.

Definition at line 223 of file MuRingForwardLayer.cc.

References theBasicComps.

Referenced by MuRingForwardDoubleLayer::subDetector().

                                                                    {
  return theBasicComps.front()->subDetector();
}

Member Data Documentation

Definition at line 65 of file MuRingForwardLayer.h.

Referenced by MuRingForwardLayer().

std::vector<const GeomDet*> MuRingForwardLayer::theBasicComps [private]

Definition at line 63 of file MuRingForwardLayer.h.

Referenced by basicComponents(), MuRingForwardLayer(), and subDetector().

Definition at line 64 of file MuRingForwardLayer.h.

Referenced by compatibleDets(), MuRingForwardLayer(), and ~MuRingForwardLayer().

Definition at line 62 of file MuRingForwardLayer.h.

Referenced by components().

std::vector<const ForwardDetRing*> MuRingForwardLayer::theRings [private]