Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:01

0001 // -*- C++ -*-
0002 //
0003 // Package:     Forward
0004 // Class  :     TotemT1Organization
0005 //
0006 // Implementation:
0007 //     <Notes on implementation>
0008 //
0009 // Original Author:  R. Capra
0010 //         Created:  Tue May 16 10:14:34 CEST 2006
0011 //
0012 
0013 // system include files
0014 
0015 // user include files
0016 #include "SimG4CMS/Forward/interface/TotemT1Organization.h"
0017 #include "SimG4CMS/Forward/interface/TotemNumberMerger.h"
0018 #include "SimG4CMS/Forward/interface/ForwardName.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 
0021 #include "G4VPhysicalVolume.hh"
0022 #include "G4VTouchable.hh"
0023 #include "globals.hh"
0024 
0025 //
0026 // constructors and destructor
0027 //
0028 TotemT1Organization ::TotemT1Organization()
0029     : _needUpdateUnitID(false),
0030       _needUpdateData(false),
0031       _currentUnitID(-1),
0032       _currentPlane(-1),
0033       _currentCSC(-1),
0034       _currentLayer(-1),
0035       _currentObjectType(Undefined) {
0036   edm::LogVerbatim("ForwardSim") << "Creating TotemT1Organization";
0037 }
0038 
0039 TotemT1Organization ::~TotemT1Organization() {}
0040 
0041 //
0042 // member functions
0043 //
0044 
0045 uint32_t TotemT1Organization ::getUnitID(const G4Step* aStep) const {
0046   return const_cast<TotemT1Organization*>(this)->getUnitID(aStep);
0047 }
0048 
0049 uint32_t TotemT1Organization ::getUnitID(const G4Step* aStep) {
0050   int currLAOT;
0051   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
0052   G4VPhysicalVolume* physVol;
0053   int ii = 0;
0054   for (ii = 0; ii < touch->GetHistoryDepth(); ii++) {
0055     physVol = touch->GetVolume(ii);
0056 
0057 #ifdef SCRIVI
0058     LogDebug("ForwardSim") << "physVol=" << physVol->GetName() << ", level=" << ii
0059                            << ", physVol->GetCopyNo()=" << physVol->GetCopyNo();
0060 #endif
0061 
0062     std::string dName = ForwardName::getName(physVol->GetName());
0063     if (dName == "TotemT1") {
0064       if (physVol->GetCopyNo() == 1)
0065         _currentDetectorPosition = 1;
0066       else if (physVol->GetCopyNo() == 2)
0067         _currentDetectorPosition = 2;
0068     }
0069   }
0070 
0071   touch = aStep->GetPreStepPoint()->GetTouchable();
0072   physVol = touch->GetVolume(0);
0073 
0074   currLAOT = physVol->GetCopyNo();
0075   _currentObjectType = static_cast<ObjectType>(currLAOT % MaxObjectTypes);
0076   _currentLayer = currLAOT / MaxObjectTypes;
0077   _currentPlane = -1;
0078   _currentCSC = -1;
0079 
0080   if (touch->GetVolume(1)) {
0081     _currentCSC = touch->GetVolume(1)->GetCopyNo();
0082     if (touch->GetVolume(2))
0083       _currentPlane = touch->GetVolume(2)->GetCopyNo();
0084   }
0085 #ifdef SCRIVI
0086   LogDebug("ForwardSim") << "CURRENT CSC " << _currentCSC << "\n"
0087                          << "CURRENT PLANE " << _currentPlane;
0088 #endif
0089   _needUpdateUnitID = true;
0090   return getCurrentUnitID();
0091 }
0092 
0093 int TotemT1Organization ::getCurrentUnitID(void) const {
0094   _checkUnitIDUpdate();
0095 #ifdef SCRIVI
0096   LogDebug("ForwardSim") << "getCurrentUnitID()=" << _currentUnitID;
0097   << endl;
0098 #endif
0099   return _currentUnitID;
0100 }
0101 
0102 void TotemT1Organization ::setCurrentUnitID(int currentUnitID) {
0103 #ifdef SCRIVI
0104   LogDebug("ForwardSim") << "_currentUnitID=" << currentUnitID;
0105 #endif
0106   _currentUnitID = currentUnitID;
0107   _needUpdateData = true;
0108 }
0109 
0110 int TotemT1Organization ::getCurrentDetectorPosition(void) const {
0111   _checkDataUpdate();
0112 #ifdef SCRIVI
0113   LogDebug("ForwardSim") << "getCurrentDetectorPosition()=" << _currentDetectorPosition;
0114 #endif
0115   return _currentDetectorPosition;
0116 }
0117 
0118 void TotemT1Organization ::setCurrentDetectorPosition(int currentDetectorPosition) {
0119 #ifdef SCRIVI
0120   LogDebug("ForwardSim") << "_currentDetectorPosition=" << currentDetectorPosition;
0121 #endif
0122   _currentDetectorPosition = currentDetectorPosition;
0123   _needUpdateUnitID = true;
0124 }
0125 
0126 int TotemT1Organization ::getCurrentPlane(void) const {
0127   _checkDataUpdate();
0128 
0129 #ifdef SCRIVI
0130   LogDebug("ForwardSim") << "getCurrentPlane()=" << _currentPlane;
0131 #endif
0132   return _currentPlane;
0133 }
0134 
0135 void TotemT1Organization ::setCurrentPlane(int currentPlane) {
0136 #ifdef SCRIVI
0137   LogDebug("ForwardSim") << "_currentPlane=" << currentPlane;
0138 #endif
0139   _currentPlane = currentPlane;
0140   _needUpdateUnitID = true;
0141 }
0142 
0143 int TotemT1Organization ::getCurrentCSC(void) const {
0144   _checkDataUpdate();
0145 #ifdef SCRIVI
0146   LogDebug("ForwardSim") << "getCurrentCSC()=" << _currentCSC;
0147 #endif
0148   return _currentCSC;
0149 }
0150 
0151 void TotemT1Organization ::setCurrentCSC(int currentCSC) {
0152 #ifdef SCRIVI
0153   LogDebug("ForwardSim") << "_currentCSC=" << currentCSC;
0154 #endif
0155   _currentCSC = currentCSC;
0156   _needUpdateUnitID = true;
0157 }
0158 
0159 int TotemT1Organization ::getCurrentLayer(void) const {
0160   _checkDataUpdate();
0161 #ifdef SCRIVI
0162   LogDebug("ForwardSim") << "getCurrentLayer()=" << _currentLayer;
0163 #endif
0164   return _currentLayer;
0165 }
0166 
0167 void TotemT1Organization ::setCurrentLayer(int currentLayer) {
0168 #ifdef SCRIVI
0169   LogDebug("ForwardSim") << "_currentLayer=" << currentLayer;
0170 #endif
0171   _currentLayer = currentLayer;
0172   _needUpdateUnitID = true;
0173 }
0174 
0175 TotemT1Organization::ObjectType TotemT1Organization ::getCurrentObjectType(void) const {
0176   _checkDataUpdate();
0177 #ifdef SCRIVI
0178   LogDebug("ForwardSim") << "getCurrentObjectType()=" << _currentObjectType;
0179 #endif
0180   return _currentObjectType;
0181 }
0182 
0183 void TotemT1Organization ::setCurrentObjectType(ObjectType currentObjectType) {
0184 #ifdef SCRIVI
0185   LogDebug("ForwardSim") << "_currentObjectType=" << currentObjectType;
0186 #endif
0187   _currentObjectType = currentObjectType;
0188   _needUpdateUnitID = true;
0189 }
0190 
0191 int TotemT1Organization ::fromObjectTypeToInt(ObjectType objectType) {
0192   int result(static_cast<int>(objectType));
0193   if (result < 0 || result >= MaxObjectTypes) {
0194     result = 0;
0195     edm::LogVerbatim("ForwardSim") << "Invalid ObjectType value (" << objectType << "). Now is \"Undefined\"";
0196   }
0197   return result;
0198 }
0199 
0200 int TotemT1Organization ::fromObjectTypeToInt(ObjectType objectType, int layer) {
0201   return fromObjectTypeToInt(objectType) + layer * MaxObjectTypes;
0202 }
0203 
0204 //
0205 // private member functions
0206 //
0207 
0208 void TotemT1Organization ::_checkUnitIDUpdate(void) const {
0209   if (_needUpdateUnitID) {
0210 #ifdef SCRIVI
0211     LogDebug("ForwardSim") << "UnitID update needed.";
0212 #endif
0213     const_cast<TotemT1Organization*>(this)->_FromDataToUnitID();
0214   } else {
0215 #ifdef SCRIVI
0216     LogDebug("ForwardSim") << "UnitID update not needed.";
0217 #endif
0218   }
0219 }
0220 
0221 void TotemT1Organization ::_checkDataUpdate(void) const {
0222   if (_needUpdateData) {
0223 #ifdef SCRIVI
0224     LogDebug("ForwardSim") << "Data update needed.";
0225 #endif
0226     const_cast<TotemT1Organization*>(this)->_FromUnitIDToData();
0227   } else {
0228 #ifdef SCRIVI
0229     LogDebug("ForwardSim") << "Data update not needed.";
0230 #endif
0231   }
0232 }
0233 
0234 void TotemT1Organization ::_FromUnitIDToData(void) {
0235   int currDP, currCSC, currOT, currPLA;
0236   unsigned long currPL, currLA;
0237 
0238   // currDP:  0..4 (5)
0239   // currPL:  0..infty
0240   // currCSC: 0..6 (7)
0241   // currLA:  0..infty
0242   // currOT:  0..MaxObjectTypes-1 (MaxObjectTypes)
0243 
0244   currDP = (_currentUnitID / 100000) % 5;  // 3;
0245   currCSC = (_currentUnitID / 5) % 7;
0246   currOT = (_currentUnitID / (5 * 7)) % MaxObjectTypes;
0247   currPLA = _currentUnitID / (5 * 7 * MaxObjectTypes);
0248 
0249   TotemNumberMerger splitter;
0250   splitter.Split(currPLA, currPL, currLA);
0251 
0252 #ifdef SCRIVI
0253   LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL << ", currCSC=" << currCSC
0254                          << ", currLA=" << currLA << ", currOT=" << currOT << ", currPLA=" << currPLA
0255                          << ", _currentUnitID=" << _currentUnitID;
0256 #endif
0257   _currentPlane = currPL - 1;
0258   _currentCSC = currCSC - 1;
0259   _currentLayer = currLA - 1;
0260   _currentObjectType = static_cast<ObjectType>(currOT);
0261 
0262   switch (currDP) {
0263     case 0:
0264       _currentDetectorPosition = 0;
0265       break;
0266     case 1:
0267       _currentDetectorPosition = 1;
0268       break;
0269     case 2:
0270       _currentDetectorPosition = 2;
0271       break;
0272     case 3:
0273       _currentDetectorPosition = 3;
0274       break;
0275     case 4:
0276       _currentDetectorPosition = 4;
0277       break;
0278   }
0279   _needUpdateData = false;
0280 }
0281 
0282 void TotemT1Organization ::_FromDataToUnitID(void) {
0283   int currDP, currPL, currCSC, currLA, currOT;
0284 #ifdef SCRIVI
0285   LogDebug("ForwardSim") << " CURRENT DETECTOR POSITION (0-3) " << _currentDetectorPosition;
0286 #endif
0287   switch (_currentDetectorPosition) {
0288     case 0:
0289       currDP = 0;
0290       break;
0291     case 1:
0292       currDP = 1;
0293       break;
0294     case 2:
0295       currDP = 2;
0296       break;
0297     case 3:
0298       currDP = 3;
0299       break;
0300     case 4:
0301       currDP = 4;
0302       break;
0303     default:
0304       _currentDetectorPosition = 0;
0305       currDP = 0;
0306       edm::LogVerbatim("ForwardSim") << "Invalid _currentDetectorPosition value (" << _currentDetectorPosition
0307                                      << "). Now is \"Undefined\"";
0308   }
0309 
0310   if (_currentPlane < -1) {
0311     edm::LogVerbatim("ForwardSim") << "Invalid _currentPlane value (" << _currentPlane << "). Now is -1";
0312     _currentPlane = -1;
0313   }
0314   currPL = _currentPlane + 1;
0315 
0316   if (_currentCSC < -1 || _currentCSC > 5) {
0317     edm::LogVerbatim("ForwardSim") << "Invalid _currentCSC value (" << _currentCSC << "). Now is -1";
0318     _currentCSC = -1;
0319   }
0320   currCSC = _currentCSC + 1;
0321 
0322   if (_currentLayer < -1) {
0323     edm::LogVerbatim("ForwardSim") << "Invalid _currentLayer value (" << _currentLayer << "). Now is -1";
0324     _currentLayer = -1;
0325   }
0326   currLA = _currentLayer + 1;
0327 
0328   currOT = fromObjectTypeToInt(_currentObjectType);
0329 
0330   // currDP:  0..2 (3)
0331   // currPL:  0..infty
0332   // currCSC: 0..6 (7)
0333   // currLA:  0..infty
0334   // currOT:  0..MaxObjectTypes-1 (MaxObjectTypes)
0335 
0336   TotemNumberMerger merger;
0337   int currPLA(merger.Merge(currPL, currLA));
0338 
0339   _currentUnitID = currDP * 100000 + 5 * (currCSC + 7 * (currOT + MaxObjectTypes * (currPLA)));
0340 #ifdef SCRIVI
0341   LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL << ", currCSC=" << currCSC
0342                          << ", currLA=" << currLA << ", currOT=" << currOT << ", currPLA=" << currPLA
0343                          << ", _currentUnitID=" << _currentUnitID;
0344 #endif
0345 
0346   _needUpdateUnitID = false;
0347 }