CMS 3D CMS Logo

RPCConeBuilder.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RPCConeBuilder
4 // Class: RPCConeBuilder
5 //
13 //
14 // Original Author: Tomasz Maciej Frueboes
15 // Created: Fri Feb 22 13:57:06 CET 2008
16 //
17 //
18 
20 
27 
28 #include <cmath>
29 #include <vector>
30 
32  m_towerBeg(iConfig.getParameter<int>("towerBeg")),
33  m_towerEnd(iConfig.getParameter<int>("towerEnd"))
34 {
35  setWhatProduced(this);
36 }
37 
38 // ------------ method called to produce the data ------------
41 {
42  auto pL1RPCConeBuilder = std::make_unique<L1RPCConeBuilder>();
43 
44  pL1RPCConeBuilder->setFirstTower(m_towerBeg);
45  pL1RPCConeBuilder->setLastTower(m_towerEnd);
46 
47  edm::ESHandle<RPCGeometry> rpcGeometry;
48  iRecord.getRecord<MuonGeometryRecord>().get(rpcGeometry);
49 
50  edm::ESHandle<L1RPCConeDefinition> l1RPCConeDefinition;
51  iRecord.getRecord<L1RPCConeDefinitionRcd>().get(l1RPCConeDefinition);
52 
54 
55  buildCones(rpcGeometry.product(), l1RPCConeDefinition.product(), ringsMap);
56 
57  // Compress all connections. Since members of this class are shared
58  // pointers this call will compress all data
59  ringsMap.begin()->second.compressConnections();
60 
61  pL1RPCConeBuilder->setConeConnectionMap(ringsMap.begin()->second.getConnectionsMap());
62 
63  pL1RPCConeBuilder->setCompressedConeConnectionMap(
64  ringsMap.begin()->second.getCompressedConnectionsMap());
65 
66  return pL1RPCConeBuilder;
67 }
68 
70  L1RPCConeDefinition const* l1RPCConeDefinition,
71  RPCStripsRing::TIdToRindMap& ringsMap) {
72 
73  // fetch geometrical data
74  auto uncompressedCons = std::make_shared<L1RPCConeBuilder::TConMap>();
75 
76  int rolls = 0;
77  for (auto const& it : rpcGeom->dets()) {
78 
79  if ( dynamic_cast< RPCRoll const * >( it ) == nullptr ) {
80  continue;
81  }
82 
83  ++rolls;
84  RPCRoll const* roll = dynamic_cast< RPCRoll const*>( it );
85 
86  int ringId = RPCStripsRing::getRingId(roll);
87  if ( ringsMap.find(ringId) == ringsMap.end() ) {
88  ringsMap[ringId] = RPCStripsRing(roll, uncompressedCons);
89  } else {
90  ringsMap[ringId].addRoll(roll);
91  }
92  }
93 
94  // filtermixed strips, fill gaps with virtual strips
95  for (auto& it : ringsMap) {
96 
97  it.second.filterOverlapingChambers();
98  it.second.fillWithVirtualStrips();
99  }
100 
101  // Xcheck, if rings are symettrical
102  for (auto& it : ringsMap) {
103  int key = it.first;
104  int sign = key/100 - (key/1000)*10;
105 
106  if (sign == 0) {
107  key += 100;
108  } else {
109  key -= 100;
110  }
111 
112  if (key != 2000){// Hey 2100 has no counterring
113  if (it.second.size() != ringsMap[key].size())
114  {
115  throw cms::Exception("RPCInternal") << " Size differs for ring " << key << " +- 100 \n";
116  }
117  }
118  }
119  buildConnections(l1RPCConeDefinition, ringsMap);
120 }
121 
123  RPCStripsRing::TIdToRindMap& ringsMap) {
124 
125  RPCStripsRing::TIdToRindMap::iterator itRef = ringsMap.begin();
126  for (;itRef != ringsMap.end(); ++itRef){ // iterate over reference rings
127 
128  RPCStripsRing::TOtherConnStructVec ringsToConnect;
129 
130  if (!itRef->second.isReferenceRing()) continue; // iterate over reference rings
131 
132  RPCStripsRing::TIdToRindMap::iterator itOther = ringsMap.begin();
133  for (;itOther != ringsMap.end(); ++itOther){ // iterate over nonreference rings
134 
135  if (itOther->second.isReferenceRing()) continue; // iterate over nonreference rings
136 
137  std::pair<int,int> pr = areConnected(itRef, itOther, l1RPCConeDefinition);
138  if ( pr.first != -1 ) {
139  RPCStripsRing::TOtherConnStruct newOtherConn;
140  newOtherConn.m_it = itOther;
141  newOtherConn.m_logplane = pr.first;
142  newOtherConn.m_logplaneSize = pr.second;
143  ringsToConnect.push_back(newOtherConn);
144  }
145  } // OtherRings iteration ends
146 
147  std::pair<int,int> prRef = areConnected(itRef, itRef, l1RPCConeDefinition);
148  if (prRef.first == -1){
149  throw cms::Exception("RPCConfig") << " Cannot determine logplane for reference ring "
150  << itRef->first << "\n ";
151  }
152 
153  itRef->second.createRefConnections(ringsToConnect, prRef.first, prRef.second);
154 
155  } // RefRings iteration ends
156 }
157 
158 // first - logplane
159 // second - logplanesize
160 std::pair<int, int>
161 RPCConeBuilder::areConnected(RPCStripsRing::TIdToRindMap::iterator ref,
162  RPCStripsRing::TIdToRindMap::iterator other,
163  L1RPCConeDefinition const* l1RPCConeDefinition) {
164 
165  int logplane = -1;
166 
167  // Do not connect rolls lying on the oposite side of detector
168  if ( ref->second.getEtaPartition()*other->second.getEtaPartition()<0 )
169  return std::make_pair(-1,0);
170 
171  int refTowerCnt = 0;
172  int index = -1;
173  int refTower = -1;
174 
175  for (auto const& itRef : l1RPCConeDefinition->getRingToTowerVec()) {
176 
177  if ( itRef.m_etaPart != std::abs(ref->second.getEtaPartition())
178  || itRef.m_hwPlane != std::abs(ref->second.getHwPlane()-1) // -1?
179  ) {
180  continue;
181  }
182 
183  ++refTowerCnt;
184  refTower = itRef.m_tower;
185 
186  for (auto const& itOther : l1RPCConeDefinition->getRingToTowerVec()) {
187 
188  if ( itOther.m_etaPart != std::abs(other->second.getEtaPartition())
189  || itOther.m_hwPlane != std::abs(other->second.getHwPlane()-1) // -1?
190  ) {
191  continue;
192  }
193 
194  if (itOther.m_tower == refTower) {
195  index = itOther.m_index;
196  }
197  }
198  }
199 
200  if (refTowerCnt > 1) {
201  throw cms::Exception("RPCConeBuilder") << " Reference(?) ring "
202  << ref->first << " "
203  << "wants to be connected to " << refTowerCnt << " towers \n";
204  }
205 
206  if (refTowerCnt == 0) {
207  throw cms::Exception("RPCConeBuilder") << " Reference(?) ring "
208  << ref->first << " "
209  << " is not connected anywhere \n";
210  }
211 
212  int lpSize = 0;
213  if (index != -1) {
214 
215  for (auto const& it : l1RPCConeDefinition->getRingToLPVec()) {
216  if (it.m_etaPart != std::abs(other->second.getEtaPartition())
217  || it.m_hwPlane != std::abs(other->second.getHwPlane()-1)
218  || it.m_index != index) {
219  continue;
220  }
221  logplane = it.m_LP;
222  }
223 
224  for (auto const& it : l1RPCConeDefinition->getLPSizeVec()) {
225  if (it.m_tower != std::abs(refTower) || it.m_LP != logplane-1) {
226  continue;
227  }
228  lpSize = it.m_size;
229  }
230 
231  //FIXME
232  if (lpSize==-1) {
233  //throw cms::Exception("getLogStrip") << " lpSize==-1\n";
234  }
235  }
236  return std::make_pair(logplane, lpSize);
237 }
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
RPCConeBuilder(const edm::ParameterSet &)
const TRingToLPVec & getRingToLPVec() const
#define nullptr
const TRingToTowerVec & getRingToTowerVec() const
std::vector< TOtherConnStruct > TOtherConnStructVec
Definition: RPCStripsRing.h:53
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const TLPSizeVec & getLPSizeVec() const
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
Definition: RPCGeometry.cc:33
std::unique_ptr< L1RPCConeBuilder > ReturnType
void buildConnections(L1RPCConeDefinition const *, RPCStripsRing::TIdToRindMap &)
std::map< int, RPCStripsRing > TIdToRindMap
Definition: RPCStripsRing.h:43
void buildCones(RPCGeometry const *, L1RPCConeDefinition const *, RPCStripsRing::TIdToRindMap &)
ReturnType produce(const L1RPCConeBuilderRcd &)
TIdToRindMap::iterator m_it
Definition: RPCStripsRing.h:50
std::pair< int, int > areConnected(RPCStripsRing::TIdToRindMap::iterator ref, RPCStripsRing::TIdToRindMap::iterator other, L1RPCConeDefinition const *)
T const * product() const
Definition: ESHandle.h:86