CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TotemT1Organization.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Forward
4 // Class : TotemT1Organization
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: R. Capra
10 // Created: Tue May 16 10:14:34 CEST 2006
11 //
12 
13 // system include files
14 
15 // user include files
19 
20 #include "G4VPhysicalVolume.hh"
21 #include "G4VTouchable.hh"
22 
23 //
24 // constructors and destructor
25 //
27  _needUpdateData(false),
28  _currentUnitID(-1),
29  _currentPlane(-1),
30  _currentCSC(-1),
31  _currentLayer(-1),
32  _currentObjectType(Undefined) {
33 
34  edm::LogInfo("ForwardSim") << "Creating TotemT1Organization";
35 }
36 
38 }
39 
40 //
41 // member functions
42 //
43 
44 uint32_t TotemT1Organization :: GetUnitID(const G4Step* aStep) const {
45  return const_cast<TotemT1Organization *>(this)->GetUnitID(aStep);
46 }
47 
48 
49 uint32_t TotemT1Organization :: GetUnitID(const G4Step* aStep) {
50 
51  int currLAOT;
52  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
53  G4VPhysicalVolume* physVol;
54  int ii =0;
55  for ( ii = 0; ii < touch->GetHistoryDepth(); ii++ ){
56  physVol = touch->GetVolume(ii);
57 
58 #ifdef SCRIVI
59  LogDebug("ForwardSim") << "physVol=" << physVol->GetName()
60  << ", level=" << ii << ", physVol->GetCopyNo()="
61  << physVol->GetCopyNo();
62 #endif
63 
64  if(physVol->GetName() == "TotemT1" &&
65  physVol->GetCopyNo()==1) _currentDetectorPosition = 1;
66  if(physVol->GetName() == "TotemT1" &&
67  physVol->GetCopyNo()==2) _currentDetectorPosition = 2;
68  }
69 
70  touch = aStep->GetPreStepPoint()->GetTouchable();
71  physVol=touch->GetVolume(0);
72 
73  currLAOT=physVol->GetCopyNo();
74  _currentObjectType=static_cast<ObjectType>(currLAOT%MaxObjectTypes);
76  _currentPlane=-1;
77  _currentCSC=-1;
78 
79  if (touch->GetVolume(1)) {
80  _currentCSC=touch->GetVolume(1)->GetCopyNo();
81  if (touch->GetVolume(2)) _currentPlane=touch->GetVolume(2)->GetCopyNo();
82  }
83 #ifdef SCRIVI
84  LogDebug("ForwardSim") << "CURRENT CSC "<<_currentCSC << "\n"
85  << "CURRENT PLANE "<<_currentPlane;
86 #endif
87  _needUpdateUnitID=true;
88  return GetCurrentUnitID();
89 }
90 
92 
94 #ifdef SCRIVI
95  LogDebug("ForwardSim") << "GetCurrentUnitID()=" << _currentUnitID;
96  << endl;
97 #endif
98  return _currentUnitID;
99 }
100 
102 
103 #ifdef SCRIVI
104  LogDebug("ForwardSim") << "_currentUnitID=" << currentUnitID;
105 #endif
106  _currentUnitID=currentUnitID;
107  _needUpdateData=true;
108 }
109 
111 
113 #ifdef SCRIVI
114  LogDebug("ForwardSim") << "GetCurrentDetectorPosition()="
116 #endif
118 }
119 
120 void TotemT1Organization :: SetCurrentDetectorPosition(int currentDetectorPosition) {
121 
122 #ifdef SCRIVI
123  LogDebug("ForwardSim") << "_currentDetectorPosition=" << currentDetectorPosition;
124 #endif
125  _currentDetectorPosition=currentDetectorPosition;
126  _needUpdateUnitID=true;
127 }
128 
130 
132 
133 #ifdef SCRIVI
134  LogDebug("ForwardSim") << "GetCurrentPlane()=" << _currentPlane;
135 #endif
136  return _currentPlane;
137 }
138 
140 
141 #ifdef SCRIVI
142  LogDebug("ForwardSim") << "_currentPlane=" << currentPlane;
143 #endif
144  _currentPlane=currentPlane;
145  _needUpdateUnitID=true;
146 }
147 
149 
151 #ifdef SCRIVI
152  LogDebug("ForwardSim") << "GetCurrentCSC()=" << _currentCSC;
153 #endif
154  return _currentCSC;
155 }
156 
158 
159 #ifdef SCRIVI
160  LogDebug("ForwardSim") << "_currentCSC=" << currentCSC;
161 #endif
162  _currentCSC=currentCSC;
163  _needUpdateUnitID=true;
164 }
165 
167 
169 #ifdef SCRIVI
170  LogDebug("ForwardSim") << "GetCurrentLayer()=" << _currentLayer;
171 #endif
172  return _currentLayer;
173 }
174 
176 
177 #ifdef SCRIVI
178  LogDebug("ForwardSim") << "_currentLayer=" << currentLayer;
179 #endif
180  _currentLayer=currentLayer;
181  _needUpdateUnitID=true;
182 }
183 
185 
187 #ifdef SCRIVI
188  LogDebug("ForwardSim") << "GetCurrentObjectType()=" << _currentObjectType;
189 #endif
190  return _currentObjectType;
191 }
192 
194 
195 #ifdef SCRIVI
196  LogDebug("ForwardSim") << "_currentObjectType=" << currentObjectType;
197 #endif
198  _currentObjectType=currentObjectType;
199  _needUpdateUnitID=true;
200 }
201 
203 
204  int result(static_cast<int>(objectType));
205  if (result<0 || result>=MaxObjectTypes) {
206  result = 0;
207  edm::LogInfo("ForwardSim") << "Invalid ObjectType value (" << objectType
208  << "). Now is \"Undefined\"";
209  }
210  return result;
211 }
212 
214  return FromObjectTypeToInt(objectType)+layer*MaxObjectTypes;
215 }
216 
217 
218 //
219 // private member functions
220 //
221 
223  if (_needUpdateUnitID) {
224 #ifdef SCRIVI
225  LogDebug("ForwardSim") << "UnitID update needed.";
226 #endif
227  const_cast<TotemT1Organization *>(this)->_FromDataToUnitID();
228  } else {
229 #ifdef SCRIVI
230  LogDebug("ForwardSim") << "UnitID update not needed.";
231 #endif
232  }
233 }
234 
236 
237  if (_needUpdateData) {
238 #ifdef SCRIVI
239  LogDebug("ForwardSim") << "Data update needed.";
240 #endif
241  const_cast<TotemT1Organization *>(this)->_FromUnitIDToData();
242  } else {
243 #ifdef SCRIVI
244  LogDebug("ForwardSim") << "Data update not needed.";
245 #endif
246  }
247 }
248 
250 
251  int currDP, currCSC, currOT, currPLA;
252  unsigned long currPL, currLA;
253 
254  // currDP: 0..4 (5)
255  // currPL: 0..infty
256  // currCSC: 0..6 (7)
257  // currLA: 0..infty
258  // currOT: 0..MaxObjectTypes-1 (MaxObjectTypes)
259 
260  currDP = (_currentUnitID/100000) % 5;// 3;
261  currCSC = (_currentUnitID/5)%7;
262  currOT = (_currentUnitID/(5*7))%MaxObjectTypes;
263  currPLA = _currentUnitID/(5*7*MaxObjectTypes);
264 
266  splitter.Split(currPLA,currPL,currLA);
267 
268 #ifdef SCRIVI
269  LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL
270  << ", currCSC=" << currCSC << ", currLA=" << currLA
271  << ", currOT=" << currOT << ", currPLA=" << currPLA
272  << ", _currentUnitID=" << _currentUnitID;
273 #endif
274  _currentPlane=currPL-1;
275  _currentCSC=currCSC-1;
276  _currentLayer=currLA-1;
277  _currentObjectType=static_cast<ObjectType>(currOT);
278 
279  switch(currDP) {
280  case 0:
282  break;
283  case 1:
285  break;
286  case 2:
288  break;
289  case 3:
291  break;
292  case 4:
294  break;
295  }
296  _needUpdateData=false;
297 }
298 
300  int currDP, currPL, currCSC, currLA, currOT;
301 #ifdef SCRIVI
302  LogDebug("ForwardSim") << " CURRENT DETECTOR POSITION (0-3) "
304 #endif
305  switch(_currentDetectorPosition) {
306  case 0:
307  currDP=0;
308  break;
309  case 1:
310  currDP=1;
311  break;
312  case 2:
313  currDP=2;
314  break;
315  case 3:
316  currDP=3;
317  break;
318  case 4:
319  currDP=4;
320  break;
321  default:
323  currDP=0;
324  edm::LogInfo("ForwardSim") << "Invalid _currentDetectorPosition value ("
326  << "). Now is \"Undefined\"";
327  }
328 
329  if (_currentPlane<-1) {
330  edm::LogInfo("ForwardSim") << "Invalid _currentPlane value ("
331  << _currentPlane << "). Now is -1";
332  _currentPlane=-1;
333  }
334  currPL=_currentPlane+1;
335 
336  if (_currentCSC<-1 || _currentCSC>5) {
337  edm::LogInfo("ForwardSim") << "Invalid _currentCSC value (" << _currentCSC
338  << "). Now is -1";
339  _currentCSC=-1;
340  }
341  currCSC=_currentCSC+1;
342 
343  if (_currentLayer<-1) {
344  edm::LogInfo("ForwardSim") << "Invalid _currentLayer value ("
345  << _currentLayer << "). Now is -1";
346  _currentLayer=-1;
347  }
348  currLA=_currentLayer+1;
349 
351 
352  // currDP: 0..2 (3)
353  // currPL: 0..infty
354  // currCSC: 0..6 (7)
355  // currLA: 0..infty
356  // currOT: 0..MaxObjectTypes-1 (MaxObjectTypes)
357 
358  TotemNumberMerger merger;
359  int currPLA(merger.Merge(currPL,currLA));
360 
361  _currentUnitID=currDP*100000+5*(currCSC+7*(currOT+MaxObjectTypes*(currPLA)));
362 #ifdef SCRIVI
363  LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL
364  << ", currCSC=" << currCSC << ", currLA=" << currLA
365  << ", currOT=" << currOT << ", currPLA=" << currPLA
366  << ", _currentUnitID=" << _currentUnitID;
367 #endif
368 
369  _needUpdateUnitID=false;
370 }
#define LogDebug(id)
void SetCurrentUnitID(int currentUnitID)
void _checkDataUpdate(void) const
int GetCurrentDetectorPosition(void) const
ObjectType GetCurrentObjectType(void) const
int GetCurrentCSC(void) const
void SetCurrentCSC(int currentCSC)
void SetCurrentLayer(int currentLayer)
unsigned long Merge(unsigned long value1, unsigned long value2) const
int ii
Definition: cuy.py:588
void Split(unsigned long source, unsigned long &value1, unsigned long &value2) const
tuple result
Definition: mps_fire.py:95
uint32_t GetUnitID(const G4Step *aStep)
def splitter
Definition: confdb.py:11
int GetCurrentUnitID(void) const
void _checkUnitIDUpdate(void) const
void SetCurrentDetectorPosition(int currentDetectorPosition)
void SetCurrentObjectType(ObjectType currentObjectType)
int GetCurrentLayer(void) const
void SetCurrentPlane(int currentPlane)
int GetCurrentPlane(void) const
volatile std::atomic< bool > shutdown_flag false
int FromObjectTypeToInt(ObjectType objectType)