![]() |
![]() |
#include <RecoLocalMuon/DTSegment/src/DTCombinatorialPatternReco4D.h>
Definition at line 42 of file DTCombinatorialPatternReco4D.h.
DTCombinatorialPatternReco4D::DTCombinatorialPatternReco4D | ( | const edm::ParameterSet & | pset | ) |
Constructor.
Definition at line 34 of file DTCombinatorialPatternReco4D.cc.
References allDTRecHits, debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), the2DAlgo, and theUpdator.
00034 : 00035 DTRecSegment4DBaseAlgo(pset), theAlgoName("DTCombinatorialPatternReco4D"){ 00036 00037 // debug parameter 00038 debug = pset.getUntrackedParameter<bool>("debug"); 00039 00040 // the updator 00041 theUpdator = new DTSegmentUpdator(pset); 00042 00043 // the input type. 00044 // If true the instructions in setDTRecSegment2DContainer will be schipped and the 00045 // theta segment will be recomputed from the 1D rechits 00046 // If false the theta segment will be taken from the Event. Caveat: in this case the 00047 // event must contain the 2D segments! 00048 allDTRecHits = pset.getParameter<bool>("AllDTRecHits"); 00049 00050 // Get the concrete 2D-segments reconstruction algo from the factory 00051 // For the 2D reco I use this reconstructor! 00052 the2DAlgo = new DTCombinatorialPatternReco(pset.getParameter<ParameterSet>("Reco2DAlgoConfig")); 00053 }
DTCombinatorialPatternReco4D::~DTCombinatorialPatternReco4D | ( | ) | [virtual] |
Destructor.
Definition at line 55 of file DTCombinatorialPatternReco4D.cc.
References the2DAlgo, and theUpdator.
00055 { 00056 delete the2DAlgo; 00057 delete theUpdator; 00058 }
virtual std::string DTCombinatorialPatternReco4D::algoName | ( | void | ) | const [inline, virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 55 of file DTCombinatorialPatternReco4D.h.
References theAlgoName.
00055 { return theAlgoName; }
std::vector<DTSegmentCand*> DTCombinatorialPatternReco4D::buildPhiSuperSegmentsCandidates | ( | std::vector< DTHitPairForFit * > & | pairPhiOwned | ) | [private] |
Referenced by reconstruct().
OwnVector< DTRecSegment4D > DTCombinatorialPatternReco4D::reconstruct | ( | void | ) | [virtual] |
Operations.
4d segment: I have the pos along the wire => further update!
Implements DTRecSegment4DBaseAlgo.
Definition at line 123 of file DTCombinatorialPatternReco4D.cc.
References allDTRecHits, edm::OwnVector< T, P >::begin(), buildPhiSuperSegmentsCandidates(), GenMuonPlsPt100GeV_cfg::cout, debug, edm::OwnVector< T, P >::end(), lat::endl(), DTChamber::id(), phi, edm::OwnVector< T, P >::push_back(), DTCombinatorialPatternReco::reconstruct(), HLT_VtxMuL3::result, sl, DTChamber::superLayer(), the2DAlgo, theChamber, theHitsFromTheta, theSegments2DTheta, theUpdator, GeomDet::toGlobal(), GeomDet::toLocal(), and DTSegmentUpdator::update().
00123 { 00124 00125 OwnVector<DTRecSegment4D> result; 00126 00127 if (debug){ 00128 cout << "Segments in " << theChamber->id() << endl; 00129 cout<<"Reconstructing of the Phi segments"<<endl; 00130 } 00131 00132 vector<DTHitPairForFit*> pairPhiOwned; 00133 vector<DTSegmentCand*> resultPhi = buildPhiSuperSegmentsCandidates(pairPhiOwned); 00134 00135 if (debug) cout << "There are " << resultPhi.size() << " Phi cand" << endl; 00136 00137 if (allDTRecHits){ 00138 // take the theta SL of this chamber 00139 const DTSuperLayer* sl = theChamber->superLayer(2); 00140 // sl points to 0 if the theta SL was not found 00141 if(sl){ 00142 // reconstruct the theta segments 00143 if(debug) cout<<"Reconstructing of the Theta segments"<<endl; 00144 OwnVector<DTSLRecSegment2D> thetaSegs = the2DAlgo->reconstruct(sl, theHitsFromTheta); 00145 vector<DTSLRecSegment2D> segments2DTheta(thetaSegs.begin(),thetaSegs.end()); 00146 theSegments2DTheta = segments2DTheta; 00147 } 00148 } 00149 00150 bool hasZed=false; 00151 00152 // has this chamber the Z-superlayer? 00153 if (theSegments2DTheta.size()){ 00154 hasZed = theSegments2DTheta.size()>0; 00155 if (debug) cout << "There are " << theSegments2DTheta.size() << " Theta cand" << endl; 00156 } else { 00157 if (debug) cout << "No Theta SL" << endl; 00158 } 00159 00160 // Now I want to build the concrete DTRecSegment4D. 00161 if(debug) cout<<"Building of the concrete DTRecSegment4D"<<endl; 00162 if (resultPhi.size()) { 00163 for (vector<DTSegmentCand*>::const_iterator phi=resultPhi.begin(); 00164 phi!=resultPhi.end(); ++phi) { 00165 00166 std::auto_ptr<DTChamberRecSegment2D> superPhi(**phi); 00167 00168 theUpdator->update(superPhi.get()); 00169 00170 if (hasZed) { 00171 00172 // Create all the 4D-segment combining the Z view with the Phi one 00173 // loop over the Z segments 00174 for(vector<DTSLRecSegment2D>::const_iterator zed = theSegments2DTheta.begin(); 00175 zed != theSegments2DTheta.end(); ++zed){ 00176 00177 // Important!! 00178 DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId()); 00179 00180 const LocalPoint posZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition() )) ; 00181 const LocalVector dirZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection() )) ; 00182 00183 DTRecSegment4D* newSeg = new DTRecSegment4D(*superPhi,*zed,posZInCh,dirZInCh); 00184 //<< 00185 00187 theUpdator->update(newSeg); 00188 if (debug) cout << "Created a 4D seg " << endl; 00189 result.push_back(newSeg); 00190 } 00191 } else { 00192 // Only phi 00193 DTRecSegment4D* newSeg = new DTRecSegment4D(*superPhi); 00194 00195 if (debug) cout << "Created a 4D segment using only the 2D Phi segment" << endl; 00196 result.push_back(newSeg); 00197 } 00198 } 00199 } else { 00200 // DTRecSegment4D from zed projection only (unlikely, not so useful, but...) 00201 if (hasZed) { 00202 for(vector<DTSLRecSegment2D>::const_iterator zed = theSegments2DTheta.begin(); 00203 zed != theSegments2DTheta.end(); ++zed){ 00204 00205 // Important!! 00206 DTSuperLayerId ZedSegSLId(zed->geographicalId().rawId()); 00207 00208 const LocalPoint posZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localPosition() )) ; 00209 const LocalVector dirZInCh = theChamber->toLocal( theChamber->superLayer(ZedSegSLId)->toGlobal(zed->localDirection() )) ; 00210 00211 DTRecSegment4D* newSeg = new DTRecSegment4D( *zed,posZInCh,dirZInCh); 00212 // << 00213 00214 if (debug) cout << "Created a 4D segment using only the 2D Theta segment" << endl; 00215 result.push_back(newSeg); 00216 } 00217 } 00218 } 00219 // finally delete the candidates! 00220 for (vector<DTSegmentCand*>::iterator phi=resultPhi.begin(); 00221 phi!=resultPhi.end(); ++phi) delete *phi; 00222 for (vector<DTHitPairForFit*>::iterator phiPair = pairPhiOwned.begin(); 00223 phiPair!=pairPhiOwned.end(); ++phiPair) delete *phiPair; 00224 return result; 00225 }
void DTCombinatorialPatternReco4D::setChamber | ( | const DTChamberId & | chId | ) | [virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 66 of file DTCombinatorialPatternReco4D.cc.
References theChamber, and theDTGeometry.
00066 { 00067 // Set the chamber 00068 theChamber = theDTGeometry->chamber(chId); 00069 }
void DTCombinatorialPatternReco4D::setDTRecHit1DContainer | ( | edm::Handle< DTRecHitCollection > | all1DHits | ) | [virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 71 of file DTCombinatorialPatternReco4D.cc.
References allDTRecHits, GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), DTChamber::id(), DTRangeMapAccessor::layersBySuperLayer(), range, theChamber, theHitsFromPhi1, theHitsFromPhi2, and theHitsFromTheta.
00071 { 00072 theHitsFromPhi1.clear(); 00073 theHitsFromPhi2.clear(); 00074 theHitsFromTheta.clear(); 00075 00076 DTRecHitCollection::range rangeHitsFromPhi1 = 00077 all1DHits->get(DTRangeMapAccessor::layersBySuperLayer( DTSuperLayerId(theChamber->id(),1) ) ); 00078 DTRecHitCollection::range rangeHitsFromPhi2 = 00079 all1DHits->get(DTRangeMapAccessor::layersBySuperLayer( DTSuperLayerId(theChamber->id(),3) ) ); 00080 00081 vector<DTRecHit1DPair> hitsFromPhi1(rangeHitsFromPhi1.first,rangeHitsFromPhi1.second); 00082 vector<DTRecHit1DPair> hitsFromPhi2(rangeHitsFromPhi2.first,rangeHitsFromPhi2.second); 00083 if(debug) 00084 cout<< "Number of DTRecHit1DPair in the SL 1 (Phi 1): " << hitsFromPhi1.size()<<endl 00085 << "Number of DTRecHit1DPair in the SL 3 (Phi 2): " << hitsFromPhi2.size()<<endl; 00086 00087 theHitsFromPhi1 = hitsFromPhi1; 00088 theHitsFromPhi2 = hitsFromPhi2; 00089 00090 if(allDTRecHits){ 00091 DTRecHitCollection::range rangeHitsFromTheta = 00092 all1DHits->get(DTRangeMapAccessor::layersBySuperLayer( DTSuperLayerId(theChamber->id(),2) ) ); 00093 00094 vector<DTRecHit1DPair> hitsFromTheta(rangeHitsFromTheta.first,rangeHitsFromTheta.second); 00095 if(debug) 00096 cout<< "Number of DTRecHit1DPair in the SL 2 (Theta): " << hitsFromTheta.size()<<endl; 00097 theHitsFromTheta = hitsFromTheta; 00098 } 00099 00100 }
void DTCombinatorialPatternReco4D::setDTRecSegment2DContainer | ( | edm::Handle< DTRecSegment2DCollection > | all2DSegments | ) | [virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 102 of file DTCombinatorialPatternReco4D.cc.
References allDTRecHits, GenMuonPlsPt100GeV_cfg::cout, debug, lat::endl(), DTChamber::id(), range, theChamber, and theSegments2DTheta.
00102 { 00103 theSegments2DTheta.clear(); 00104 00105 if(!allDTRecHits){ 00106 00107 //Extract the DTRecSegment2DCollection range for the theta SL 00108 DTRecSegment2DCollection::range rangeTheta = 00109 all2DSegments->get(DTSuperLayerId(theChamber->id(),2)); 00110 00111 // Fill the DTRecSegment2D container for the theta SL 00112 vector<DTSLRecSegment2D> segments2DTheta(rangeTheta.first,rangeTheta.second); 00113 00114 if(debug) 00115 cout << "Number of 2D-segments in the second SL (Theta): " << segments2DTheta.size() << endl; 00116 theSegments2DTheta = segments2DTheta; 00117 } 00118 00119 }
void DTCombinatorialPatternReco4D::setES | ( | const edm::EventSetup & | setup | ) | [virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 60 of file DTCombinatorialPatternReco4D.cc.
References edm::EventSetup::get(), DTCombinatorialPatternReco::setES(), DTSegmentUpdator::setES(), the2DAlgo, theDTGeometry, and theUpdator.
00060 { 00061 setup.get<MuonGeometryRecord>().get(theDTGeometry); 00062 the2DAlgo->setES(setup); 00063 theUpdator->setES(setup); 00064 }
virtual bool DTCombinatorialPatternReco4D::wants2DSegments | ( | ) | [inline, virtual] |
Implements DTRecSegment4DBaseAlgo.
Definition at line 61 of file DTCombinatorialPatternReco4D.h.
References allDTRecHits.
00061 {return !allDTRecHits;}
Definition at line 86 of file DTCombinatorialPatternReco4D.h.
Referenced by DTCombinatorialPatternReco4D(), reconstruct(), setDTRecHit1DContainer(), setDTRecSegment2DContainer(), and wants2DSegments().
bool DTCombinatorialPatternReco4D::debug [private] |
Definition at line 70 of file DTCombinatorialPatternReco4D.h.
Referenced by DTCombinatorialPatternReco4D(), reconstruct(), setDTRecHit1DContainer(), and setDTRecSegment2DContainer().
Definition at line 78 of file DTCombinatorialPatternReco4D.h.
Referenced by DTCombinatorialPatternReco4D(), reconstruct(), setES(), and ~DTCombinatorialPatternReco4D().
std::string DTCombinatorialPatternReco4D::theAlgoName [private] |
const DTChamber* DTCombinatorialPatternReco4D::theChamber [private] |
Definition at line 83 of file DTCombinatorialPatternReco4D.h.
Referenced by reconstruct(), setChamber(), setDTRecHit1DContainer(), and setDTRecSegment2DContainer().
Definition at line 74 of file DTCombinatorialPatternReco4D.h.
Referenced by setChamber(), and setES().
std::vector<DTRecHit1DPair> DTCombinatorialPatternReco4D::theHitsFromPhi1 [private] |
Definition at line 90 of file DTCombinatorialPatternReco4D.h.
Referenced by setDTRecHit1DContainer().
std::vector<DTRecHit1DPair> DTCombinatorialPatternReco4D::theHitsFromPhi2 [private] |
Definition at line 92 of file DTCombinatorialPatternReco4D.h.
Referenced by setDTRecHit1DContainer().
std::vector<DTRecHit1DPair> DTCombinatorialPatternReco4D::theHitsFromTheta [private] |
Definition at line 91 of file DTCombinatorialPatternReco4D.h.
Referenced by reconstruct(), and setDTRecHit1DContainer().
std::vector<DTSLRecSegment2D> DTCombinatorialPatternReco4D::theSegments2DTheta [private] |
Definition at line 89 of file DTCombinatorialPatternReco4D.h.
Referenced by reconstruct(), and setDTRecSegment2DContainer().
Definition at line 81 of file DTCombinatorialPatternReco4D.h.
Referenced by DTCombinatorialPatternReco4D(), reconstruct(), setES(), and ~DTCombinatorialPatternReco4D().