CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/SimG4CMS/Forward/src/TotemT1Organization.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     Forward
00004 // Class  :     TotemT1Organization
00005 //
00006 // Implementation:
00007 //     <Notes on implementation>
00008 //
00009 // Original Author:  R. Capra
00010 //         Created:  Tue May 16 10:14:34 CEST 2006
00011 // $Id: TotemT1Organization.cc,v 1.1 2006/05/17 16:18:58 sunanda Exp $
00012 //
00013 
00014 // system include files
00015 
00016 // user include files
00017 #include "SimG4CMS/Forward/interface/TotemT1Organization.h"
00018 #include "SimG4CMS/Forward/interface/TotemNumberMerger.h"
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 
00021 #include "G4VPhysicalVolume.hh"
00022 #include "G4VTouchable.hh" 
00023 
00024 //
00025 // constructors and destructor
00026 //
00027 TotemT1Organization :: TotemT1Organization() :  _needUpdateUnitID(false),
00028                                                 _needUpdateData(false),
00029                                                 _currentUnitID(-1),
00030                                                 _currentPlane(-1),
00031                                                 _currentCSC(-1),
00032                                                 _currentLayer(-1),
00033                                                 _currentObjectType(Undefined) {
00034 
00035   edm::LogInfo("ForwardSim") << "Creating TotemT1Organization";
00036 }
00037 
00038 TotemT1Organization :: ~TotemT1Organization() {
00039 }
00040 
00041 //
00042 // member functions
00043 //
00044 
00045 uint32_t TotemT1Organization :: GetUnitID(const G4Step* aStep) const {
00046   return const_cast<TotemT1Organization *>(this)->GetUnitID(aStep);
00047 }
00048 
00049 
00050 uint32_t TotemT1Organization :: GetUnitID(const G4Step* aStep) {
00051 
00052   int currLAOT;
00053   const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00054   G4VPhysicalVolume* physVol;
00055   int ii =0;
00056   for ( ii = 0; ii < touch->GetHistoryDepth(); ii++ ){
00057     physVol = touch->GetVolume(ii);
00058 
00059 #ifdef SCRIVI
00060     LogDebug("ForwardSim") << "physVol=" << physVol->GetName()
00061                            << ", level=" << ii  << ", physVol->GetCopyNo()=" 
00062                            << physVol->GetCopyNo();
00063 #endif
00064 
00065     if(physVol->GetName() == "TotemT1" && 
00066        physVol->GetCopyNo()==1) _currentDetectorPosition = 1;
00067     if(physVol->GetName() == "TotemT1" && 
00068        physVol->GetCopyNo()==2) _currentDetectorPosition = 2;
00069   }
00070 
00071   touch = aStep->GetPreStepPoint()->GetTouchable();
00072   physVol=touch->GetVolume(0);
00073  
00074   currLAOT=physVol->GetCopyNo();
00075   _currentObjectType=static_cast<ObjectType>(currLAOT%MaxObjectTypes);
00076   _currentLayer=currLAOT/MaxObjectTypes;
00077   _currentPlane=-1;
00078   _currentCSC=-1;
00079 
00080   if (touch->GetVolume(1)) {
00081     _currentCSC=touch->GetVolume(1)->GetCopyNo();
00082     if (touch->GetVolume(2)) _currentPlane=touch->GetVolume(2)->GetCopyNo();
00083   }
00084 #ifdef SCRIVI
00085   LogDebug("ForwardSim") << "CURRENT CSC "<<_currentCSC << "\n"
00086                          << "CURRENT PLANE "<<_currentPlane;
00087 #endif
00088   _needUpdateUnitID=true;
00089   return GetCurrentUnitID();
00090 }
00091 
00092 int TotemT1Organization :: GetCurrentUnitID(void) const {
00093 
00094   _checkUnitIDUpdate();
00095 #ifdef SCRIVI
00096  LogDebug("ForwardSim") << "GetCurrentUnitID()=" << _currentUnitID;
00097                                              << endl;
00098 #endif
00099  return _currentUnitID;
00100 }
00101 
00102 void TotemT1Organization :: SetCurrentUnitID(int currentUnitID) {
00103 
00104 #ifdef SCRIVI
00105   LogDebug("ForwardSim") << "_currentUnitID=" << currentUnitID;
00106 #endif
00107   _currentUnitID=currentUnitID;
00108   _needUpdateData=true;
00109 }
00110 
00111 int  TotemT1Organization :: GetCurrentDetectorPosition(void) const {
00112 
00113   _checkDataUpdate();
00114 #ifdef SCRIVI
00115   LogDebug("ForwardSim") << "GetCurrentDetectorPosition()=" 
00116                          << _currentDetectorPosition;
00117 #endif
00118  return _currentDetectorPosition;
00119 }
00120 
00121 void TotemT1Organization :: SetCurrentDetectorPosition(int currentDetectorPosition) {
00122 
00123 #ifdef SCRIVI
00124   LogDebug("ForwardSim") << "_currentDetectorPosition=" << currentDetectorPosition;
00125 #endif
00126   _currentDetectorPosition=currentDetectorPosition;
00127   _needUpdateUnitID=true;
00128 }
00129 
00130 int TotemT1Organization :: GetCurrentPlane(void) const {
00131 
00132   _checkDataUpdate();
00133 
00134 #ifdef SCRIVI
00135   LogDebug("ForwardSim") << "GetCurrentPlane()=" << _currentPlane;
00136 #endif
00137  return _currentPlane;
00138 }
00139 
00140 void TotemT1Organization :: SetCurrentPlane(int currentPlane) {
00141 
00142 #ifdef SCRIVI
00143   LogDebug("ForwardSim") << "_currentPlane=" << currentPlane;
00144 #endif
00145   _currentPlane=currentPlane;
00146   _needUpdateUnitID=true;
00147 }
00148 
00149 int TotemT1Organization :: GetCurrentCSC(void) const {
00150 
00151   _checkDataUpdate();
00152 #ifdef SCRIVI
00153   LogDebug("ForwardSim") << "GetCurrentCSC()=" << _currentCSC;
00154 #endif 
00155  return _currentCSC;
00156 }
00157 
00158 void TotemT1Organization :: SetCurrentCSC(int currentCSC) {
00159 
00160 #ifdef SCRIVI
00161   LogDebug("ForwardSim") << "_currentCSC=" << currentCSC;
00162 #endif
00163   _currentCSC=currentCSC;
00164   _needUpdateUnitID=true; 
00165 }
00166 
00167 int TotemT1Organization :: GetCurrentLayer(void) const {
00168 
00169   _checkDataUpdate();
00170 #ifdef SCRIVI
00171   LogDebug("ForwardSim") << "GetCurrentLayer()=" << _currentLayer;
00172 #endif
00173   return _currentLayer;
00174 }
00175 
00176 void TotemT1Organization :: SetCurrentLayer(int currentLayer) {
00177 
00178 #ifdef SCRIVI
00179   LogDebug("ForwardSim") << "_currentLayer=" << currentLayer;
00180 #endif
00181   _currentLayer=currentLayer;
00182   _needUpdateUnitID=true;
00183 }
00184 
00185 TotemT1Organization::ObjectType  TotemT1Organization :: GetCurrentObjectType(void) const {
00186   
00187  _checkDataUpdate();
00188 #ifdef SCRIVI
00189  LogDebug("ForwardSim") << "GetCurrentObjectType()=" << _currentObjectType;
00190 #endif
00191  return _currentObjectType;
00192 }
00193 
00194 void TotemT1Organization :: SetCurrentObjectType(ObjectType currentObjectType){
00195 
00196 #ifdef SCRIVI
00197   LogDebug("ForwardSim") << "_currentObjectType=" << currentObjectType;
00198 #endif
00199   _currentObjectType=currentObjectType;
00200   _needUpdateUnitID=true;
00201 }
00202 
00203 int TotemT1Organization :: FromObjectTypeToInt(ObjectType objectType) {
00204 
00205   int result(static_cast<int>(objectType));
00206   if (result<0 || result>=MaxObjectTypes) {
00207     result = 0;
00208     edm::LogInfo("ForwardSim") << "Invalid ObjectType value (" << objectType
00209                                << "). Now is \"Undefined\"";
00210   }
00211   return result;
00212 }
00213 
00214 int TotemT1Organization :: FromObjectTypeToInt(ObjectType objectType, int layer) {
00215   return FromObjectTypeToInt(objectType)+layer*MaxObjectTypes;
00216 }
00217 
00218 
00219 //
00220 // private member functions
00221 //
00222 
00223 void TotemT1Organization :: _checkUnitIDUpdate(void) const {
00224   if (_needUpdateUnitID) {
00225 #ifdef SCRIVI
00226     LogDebug("ForwardSim") << "UnitID update needed.";
00227 #endif
00228     const_cast<TotemT1Organization *>(this)->_FromDataToUnitID();
00229   } else {
00230 #ifdef SCRIVI
00231     LogDebug("ForwardSim") << "UnitID update not needed.";
00232 #endif
00233   }
00234 }
00235 
00236 void TotemT1Organization :: _checkDataUpdate(void) const {
00237 
00238   if (_needUpdateData) {
00239 #ifdef SCRIVI
00240     LogDebug("ForwardSim") << "Data update needed.";
00241 #endif
00242     const_cast<TotemT1Organization *>(this)->_FromUnitIDToData();
00243   } else {
00244 #ifdef SCRIVI
00245     LogDebug("ForwardSim") << "Data update not needed.";
00246 #endif
00247   }
00248 }
00249 
00250 void TotemT1Organization :: _FromUnitIDToData(void) {
00251 
00252   int currDP, currCSC, currOT, currPLA;
00253   unsigned long currPL, currLA;
00254  
00255   // currDP:  0..4 (5)
00256   // currPL:  0..infty
00257   // currCSC: 0..6 (7)
00258   // currLA:  0..infty
00259   // currOT:  0..MaxObjectTypes-1 (MaxObjectTypes)
00260 
00261   currDP  = (_currentUnitID/100000) % 5;// 3;
00262   currCSC = (_currentUnitID/5)%7;
00263   currOT  = (_currentUnitID/(5*7))%MaxObjectTypes;
00264   currPLA =  _currentUnitID/(5*7*MaxObjectTypes);
00265  
00266   TotemNumberMerger splitter;
00267   splitter.Split(currPLA,currPL,currLA);
00268  
00269 #ifdef SCRIVI
00270   LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL  
00271                          << ", currCSC=" << currCSC << ", currLA=" << currLA  
00272                          << ", currOT=" << currOT << ", currPLA=" << currPLA  
00273                          << ", _currentUnitID=" << _currentUnitID;
00274 #endif
00275   _currentPlane=currPL-1;
00276   _currentCSC=currCSC-1;
00277   _currentLayer=currLA-1;
00278   _currentObjectType=static_cast<ObjectType>(currOT);
00279 
00280   switch(currDP) {
00281   case 0: 
00282     _currentDetectorPosition=0;
00283     break;
00284   case 1: 
00285     _currentDetectorPosition=1;
00286     break;
00287   case 2: 
00288     _currentDetectorPosition=2;
00289     break;
00290   case 3:
00291     _currentDetectorPosition=3;
00292     break;
00293   case 4:
00294     _currentDetectorPosition=4;
00295     break;
00296   } 
00297   _needUpdateData=false;
00298 }
00299 
00300 void TotemT1Organization :: _FromDataToUnitID(void) {
00301   int currDP, currPL, currCSC, currLA, currOT;
00302 #ifdef SCRIVI
00303   LogDebug("ForwardSim") << " CURRENT DETECTOR POSITION (0-3) " 
00304                          <<_currentDetectorPosition;
00305 #endif
00306   switch(_currentDetectorPosition)  {
00307   case 0:
00308     currDP=0;
00309     break;
00310   case 1:
00311     currDP=1;
00312     break;
00313   case 2:
00314     currDP=2;
00315     break;
00316   case 3:
00317     currDP=3;
00318     break;
00319   case 4:
00320     currDP=4;
00321     break;
00322   default:
00323     _currentDetectorPosition=0;
00324     currDP=0;
00325     edm::LogInfo("ForwardSim") << "Invalid _currentDetectorPosition value (" 
00326                                << _currentDetectorPosition
00327                                << "). Now is \"Undefined\"";
00328   }
00329  
00330   if (_currentPlane<-1) {
00331     edm::LogInfo("ForwardSim") << "Invalid _currentPlane value (" 
00332                                << _currentPlane << "). Now is -1";
00333     _currentPlane=-1;
00334   }
00335   currPL=_currentPlane+1;
00336   
00337   if (_currentCSC<-1 || _currentCSC>5)  {
00338     edm::LogInfo("ForwardSim") << "Invalid _currentCSC value (" << _currentCSC
00339                                << "). Now is -1";
00340     _currentCSC=-1;
00341   }
00342   currCSC=_currentCSC+1;
00343   
00344   if (_currentLayer<-1)  {
00345     edm::LogInfo("ForwardSim") << "Invalid _currentLayer value (" 
00346                                << _currentLayer << "). Now is -1";
00347     _currentLayer=-1;
00348   }
00349   currLA=_currentLayer+1;
00350  
00351   currOT=FromObjectTypeToInt(_currentObjectType);
00352  
00353   // currDP:  0..2 (3)
00354   // currPL:  0..infty
00355   // currCSC: 0..6 (7)
00356   // currLA:  0..infty
00357   // currOT:  0..MaxObjectTypes-1 (MaxObjectTypes)
00358 
00359   TotemNumberMerger merger;
00360   int currPLA(merger.Merge(currPL,currLA));
00361  
00362   _currentUnitID=currDP*100000+5*(currCSC+7*(currOT+MaxObjectTypes*(currPLA)));
00363 #ifdef SCRIVI
00364   LogDebug("ForwardSim") << "currDP=" << currDP << ", currPL=" << currPL  
00365                          << ", currCSC=" << currCSC << ", currLA=" << currLA  
00366                          << ", currOT=" << currOT << ", currPLA=" << currPLA  
00367                          << ", _currentUnitID=" << _currentUnitID;
00368 #endif
00369  
00370   _needUpdateUnitID=false;
00371 }