| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

In This Package:

DsPmtSensDet Class Reference

#include <DsPmtSensDet.h>

Inheritance diagram for DsPmtSensDet:

[legend]
Collaboration diagram for DsPmtSensDet:
[legend]
List of all members.

Public Types

 SUCCESS
 NO_INTERFACE
 VERSMISMATCH
 LAST_ERROR
enum  Status
enum  Status
enum  Status

Public Member Functions

 DsPmtSensDet (const std::string &type, const std::string &name, const IInterface *parent)
virtual ~DsPmtSensDet ()
virtual void Initialize (G4HCofThisEvent *HCE)
virtual void EndOfEvent (G4HCofThisEvent *HCE)
virtual bool ProcessHits (G4Step *step, G4TouchableHistory *history)
virtual StatusCode initialize ()
virtual StatusCode finalize ()
virtual bool processStep (G4Step *step, G4TouchableHistory *history)
virtual unsigned long release ()
StatusCode release (const IInterface *interface) const
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvInterface)=0
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvInterface)=0
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvInterface)=0
virtual StatusCode queryInterface (const InterfaceID &riid, void **ppvUnknown)
virtual unsigned long addRef ()=0
virtual unsigned long addRef ()=0
virtual unsigned long addRef ()=0
virtual unsigned long addRef ()
virtual const std::string & type () const =0
virtual const std::string & type () const
virtual const IInterfaceparent () const =0
virtual const IInterfaceparent () const
virtual StatusCode configure ()=0
virtual StatusCode configure ()
virtual StatusCode start ()=0
virtual StatusCode start ()
virtual StatusCode stop ()=0
virtual StatusCode stop ()
virtual StatusCode terminate ()=0
virtual StatusCode terminate ()
virtual StatusCode reinitialize ()=0
virtual StatusCode reinitialize ()
virtual StatusCode restart ()=0
virtual StatusCode restart ()
virtual Gaudi::StateMachine::State FSMState () const =0
virtual Gaudi::StateMachine::State FSMState () const
virtual StatusCode sysInitialize ()=0
virtual StatusCode sysInitialize ()
virtual StatusCode sysStart ()=0
virtual StatusCode sysStart ()
virtual StatusCode sysStop ()=0
virtual StatusCode sysStop ()
virtual StatusCode sysFinalize ()=0
virtual StatusCode sysFinalize ()
virtual StatusCode sysReinitialize ()=0
virtual StatusCode sysReinitialize ()
virtual StatusCode sysRestart ()=0
virtual StatusCode sysRestart ()
virtual unsigned long refCount () const =0
virtual const std::string & name () const =0
virtual const std::string & name () const
virtual void handle (const Incident &i)
IGiGaSvcgigaSvc () const
IGiGaSetUpSvcsetupSvc () const
INTupleSvcntupleSvc () const
INTupleSvcevtColSvc () const
IDataProviderSvcdetSvc () const
IDataProviderSvcevtSvc () const
IIncidentSvcincSvc () const
IChronoStatSvcchronoSvc () const
IHistogramSvchistoSvc () const
IAlgContextSvccontextSvc () const
DataObjectput (IDataProviderSvc *svc, DataObject *object, const std::string &address, const bool useRootInTES=true) const
DataObjectput (DataObject *object, const std::string &address, const bool useRootInTES=true) const
Gaudi::Utils::GetData< TYPE
>::return_type 
get (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
Gaudi::Utils::GetData< TYPE
>::return_type 
get (const std::string &location, const bool useRootInTES=true) const
TYPE * getDet (IDataProviderSvc *svc, const std::string &location) const
TYPE * getDet (const std::string &location) const
bool exist (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
bool exist (const std::string &location, const bool useRootInTES=true) const
bool existDet (IDataProviderSvc *svc, const std::string &location) const
bool existDet (const std::string &location) const
TYPE * getOrCreate (IDataProviderSvc *svc, const std::string &location, const bool useRootInTES=true) const
TYPE * getOrCreate (const std::string &location, const bool useRootInTES=true) const
TOOL * tool (const std::string &type, const std::string &name, const IInterface *parent=0, bool create=true) const
TOOL * tool (const std::string &type, const IInterface *parent=0, bool create=true) const
SERVICE * svc (const std::string &name, const bool create=true) const
IUpdateManagerSvcupdMgrSvc () const
IDataProviderSvcfastContainersSvc () const
StatusCode Error (const std::string &msg, const StatusCode st=StatusCode::FAILURE, const size_t mx=10) const
StatusCode Warning (const std::string &msg, const StatusCode st=StatusCode::FAILURE, const size_t mx=10) const
StatusCode Print (const std::string &msg, const StatusCode st=StatusCode::SUCCESS, const MSG::Level lev=MSG::INFO) const
StatusCode Assert (const bool ok, const std::string &message="", const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Assert (const bool ok, const char *message, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg, const GaudiException &exc, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg, const std::exception &exc, const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
StatusCode Exception (const std::string &msg="no message", const StatusCode sc=StatusCode(StatusCode::FAILURE, true)) const
MsgStreammsgStream (const MSG::Level level) const
MsgStreamalways () const
MsgStreamfatal () const
MsgStreamerr () const
MsgStreamerror () const
MsgStreamwarning () const
MsgStreaminfo () const
MsgStreamdebug () const
MsgStreamverbose () const
MsgStreammsg () const
const Statisticscounters () const
StatEntitycounter (const std::string &tag) const
MSG::Level msgLevel () const
bool msgLevel (const MSG::Level level) const
void resetMsgStream () const
bool typePrint () const
bool propsPrint () const
bool statPrint () const
bool errorsPrint () const
long printStat (const MSG::Level level=MSG::ALWAYS) const
long printErrors (const MSG::Level level=MSG::ALWAYS) const
long printProps (const MSG::Level level=MSG::ALWAYS) const
void registerCondition (const std::string &condition, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (const std::string &condition, CondType *&condPtrDest, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (char *condition, StatusCode(CallerClass::*mf)()=NULL)
void registerCondition (TargetClass *condition, StatusCode(CallerClass::*mf)()=NULL)
StatusCode runUpdate ()
TransientFastContainer< T > * getFastContainer (const std::string &location, typename TransientFastContainer< T >::size_type initial=0)
const std::string & context () const
const std::string & rootInTES () const
double globalTimeOffset () const
virtual Gaudi::StateMachine::State targetFSMState () const
virtual StatusCode setProperty (const Property &p)
virtual StatusCode setProperty (const std::string &s)
virtual StatusCode setProperty (const std::string &n, const std::string &v)
StatusCode setProperty (const std::string &name, const TYPE &value)
virtual StatusCode getProperty (Property *p) const
virtual const PropertygetProperty (const std::string &name) const
virtual StatusCode getProperty (const std::string &n, std::string &v) const
virtual const std::vector<
Property * > & 
getProperties () const
PropertyMgrgetPropertyMgr ()
ISvcLocatorserviceLocator () const
ISvcLocatorsvcLoc () const
IMessageSvcmsgSvc () const
IToolSvctoolSvc () const
StatusCode setProperties ()
StatusCode service (const std::string &name, T *&svc, bool createIf=true) const
StatusCode service (const std::string &type, const std::string &name, T *&svc) const
void declInterface (const InterfaceID &, void *)
PropertydeclareProperty (const std::string &name, T &property, const std::string &doc="none") const
PropertydeclareRemoteProperty (const std::string &name, IProperty *rsvc, const std::string &rname="") const
IAuditorSvcauditorSvc () const
IMonitorSvcmonitorSvc () const
void declareInfo (const std::string &name, const T &var, const std::string &desc) const
void declareInfo (const std::string &name, const std::string &format, const void *var, int size, const std::string &desc) const

Static Public Member Functions

static const InterfaceIDinterfaceID ()
static const InterfaceIDinterfaceID ()
static const InterfaceIDinterfaceID ()

Public Attributes

 SUCCESS
 NO_INTERFACE
 VERSMISMATCH
 LAST_ERROR

Protected Types

typedef std::map< std::string,
StatEntity
Statistics
typedef std::map< std::string,
unsigned int > 
Counter
typedef std::vector< IAlgTool * > AlgTools
typedef std::pair< IInterface *,
std::string > 
ServiceEntry
typedef std::vector< ServiceEntryServices

Protected Member Functions

virtual unsigned long refCount () const
StatusCode releaseTool (const IAlgTool *tool) const
StatusCode releaseSvc (const IInterface *svc) const
int outputLevel () const
IntegerPropertyoutputLevelProperty ()
void initOutputLevel (Property &prop)

Static Protected Attributes

static const bool IgnoreRootInTES
static const bool UseRootInTES

Private Types

typedef std::map< short int,
G4DhHitCollection * > 
LocalHitCache

Private Member Functions

void StoreHit (DayaBay::SimPmtHit *hit, int trackid)
const DetectorElementSensDetElem (const G4TouchableHistory &hist)
int SensDetId (const DetectorElement &de)
double SensDetQE (G4LogicalVolume *logvol, double energy)

Private Attributes

std::vector< std::string > m_cathodeLogVols
 CathodeLogicalVolumes : name of logical volumes in which this sensitive detector is operating.
std::vector< std::string > m_sensorStructures
 SensorStructures : names of paths in TDS in which to search for sensor detector elements using this sensitive detector.
std::string m_idParameter
 PackedIdParameterName : name of user paramater of the counted detector element which holds the packed, globally unique PMT ID.
std::string m_t2deName
 TouchableToDetelem : the ITouchableToDetectorElement to use to resolve sensor ID.
ITouchableToDetectorElementm_t2de
double m_qeScale
 QEScale: Upward adjustment of DetSim efficiency to allow PMT-to-PMT efficiency variation in the electronics simulation.
bool m_applyPreQE
 preQE is maximal possible PMT QE.
double m_preQE
std::string m_qeffParamName
 QEffParameterName : name of user parameter in the photo cathode volume that holds the quantum efficiency tabproperty.
bool m_preserveWeight
 PreserveWeight : If true, any track weight will simply be passed on in the Hit's weight data member.
LocalHitCache m_hc
Rndm::Numbers m_uni

Detailed Description

Definition at line 26 of file DsPmtSensDet.h.


Member Typedef Documentation

typedef std::map<short int,G4DhHitCollection*> DsPmtSensDet::LocalHitCache [private]

Definition at line 95 of file DsPmtSensDet.h.


Constructor & Destructor Documentation

DsPmtSensDet::DsPmtSensDet ( const std::string &  type,
const std::string &  name,
const IInterface parent 
)

Definition at line 52 of file DsPmtSensDet.cc.

00055     : G4VSensitiveDetector(name)
00056     , GiGaSensDetBase(type,name,parent)
00057     , m_t2de(0)
00058 {
00059     info() << "DsPmtSensDet (" << type << "/" << name << ") created" << endreq;
00060 
00061     declareProperty("CathodeLogicalVolume",
00062                     m_cathodeLogVols,
00063                     "Photo-Cathode logical volume to which this SD is attached.");
00064 
00065     declareProperty("TouchableToDetelem", m_t2deName = "TH2DE",
00066                     "The ITouchableToDetectorElement to use to resolve sensor.");
00067 
00068     declareProperty("SensorStructures",m_sensorStructures,
00069                     "TDS Paths in which to look for sensor detector elements"
00070                     " using this sensitive detector");
00071 
00072     declareProperty("PackedIdPropertyName",m_idParameter="PmtID",
00073                     "The name of the user property holding the PMT ID.");
00074 
00075     declareProperty("QEffParameterName",m_qeffParamName="EFFICIENCY",
00076                     "name of user parameter in the photo cathode volume that"
00077                     " holds the quantum efficiency tabproperty");
00078 
00079     declareProperty("QEScale",m_qeScale=1.0 / 0.9,
00080                     "Upward scaling of the quantum efficiency by inverse of mean PMT-to-PMT efficiency in electronics simulation.");
00081     
00082     declareProperty("PreQE", m_preQE=0.32,
00083                     "Miximal possible PMT QE value to be applied on the stage of "
00084                     "generating photons in scintillation/Cerenkov processes");
00085 
00086     declareProperty("ApplyPreQE",m_applyPreQE=true,
00087                     "Apply preQE on stage of generating photons of not.");
00088 
00089     m_cathodeLogVols.push_back("/dd/Geometry/PMT/lvPmtHemiCathode");
00090     m_cathodeLogVols.push_back("/dd/Geometry/PMT/lvHeadonPmtCathode");
00091 }

DsPmtSensDet::~DsPmtSensDet (  )  [virtual]

Definition at line 93 of file DsPmtSensDet.cc.

00094 {
00095 }


Member Function Documentation

void DsPmtSensDet::Initialize ( G4HCofThisEvent *  HCE  )  [virtual]

Definition at line 198 of file DsPmtSensDet.cc.

00199 {
00200     m_hc.clear();
00201 
00202     G4DhHitCollection* hc = new G4DhHitCollection(SensitiveDetectorName,collectionName[0]);
00203     m_hc[0] = hc;
00204     int hcid = G4SDManager::GetSDMpointer()->GetCollectionID(hc);
00205     hce->AddHitsCollection(hcid,hc);
00206 
00207     for (int isite=0; site_ids[isite] >= 0; ++isite) {
00208         for (int idet=0; detector_ids[idet] >= 0; ++idet) {
00209             DayaBay::Detector det(site_ids[isite],detector_ids[idet]);
00210 
00211             if (det.bogus()) continue;
00212 
00213             string name=det.detName();
00214             G4DhHitCollection* hc = new G4DhHitCollection(SensitiveDetectorName,name.c_str());
00215             short int id = det.siteDetPackedData();
00216             m_hc[id] = hc;
00217 
00218             int hcid = G4SDManager::GetSDMpointer()->GetCollectionID(hc);
00219             hce->AddHitsCollection(hcid,hc);
00220             debug() << "Add hit collection with hcid=" << hcid << ", cached ID=" 
00221                     << (void*)id 
00222                     << " name= \"" << SensitiveDetectorName << "/" << name << "\"" 
00223                     << endreq;
00224         }
00225     }
00226 
00227     debug() << "DsPmtSensDet Initialize, made "
00228            << hce->GetNumberOfCollections() << " collections"
00229            << endreq;
00230     
00231 }

void DsPmtSensDet::EndOfEvent ( G4HCofThisEvent *  HCE  )  [virtual]

Definition at line 497 of file DsPmtSensDet.cc.

00498 {
00499     debug() << "Cache has " << m_hc.size() << " collections" << endreq;
00500 
00501     int ncols = hce->GetNumberOfCollections();
00502     debug() << "DsPmtSensDet EndOfEvent " << ncols << " collections.";
00503     //    for (int ind=0; ind<ncols; ++ind) {
00504     //  G4VHitsCollection* hc = hce->GetHC(ind);
00505     //  info () << ind << ": " 
00506     //          << hc->GetSDname() << "//" << hc->GetName() << " has " 
00507     //          << hc->GetSize() << " hits\n";
00508     //}
00509     //info() << endreq;
00510 
00511 
00512     int tothits = 0;
00513     for (int ind=0; ind<ncols; ++ind) {
00514       G4VHitsCollection* hc = hce->GetHC(ind);
00515       if ( hc->GetSize() > 0) 
00516         {
00517           if ( tothits == 0) debug() << endreq;
00518           debug() << ind << ": " 
00519                   << hc->GetSDname() << "//" << hc->GetName() << " has " 
00520                   << hc->GetSize() << " hits" << endreq;
00521         }
00522       tothits += hc->GetSize() ;
00523     }
00524     if ( tothits == 0 ) debug() << " No hits found in " << ncols << " collections."  << endreq;
00525 }

bool DsPmtSensDet::ProcessHits ( G4Step *  step,
G4TouchableHistory *  history 
) [virtual]

Definition at line 321 of file DsPmtSensDet.cc.

00323 {
00324     //if (!step) return false; just crash for now if not defined
00325 
00326     // Find out what detector we are in (ADx, IWS or OWS)
00327     G4StepPoint* preStepPoint = step->GetPreStepPoint();
00328 
00329     double energyDep = step->GetTotalEnergyDeposit();
00330 
00331     if (energyDep <= 0.0) {
00332         //debug() << "Hit energy too low: " << energyDep/CLHEP::eV << endreq;
00333         return false;
00334     }
00335 
00336     const G4TouchableHistory* hist = 
00337         dynamic_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
00338     if (!hist or !hist->GetHistoryDepth()) {
00339         error() << "ProcessHits: step has no or empty touchable history" << endreq;
00340         return false;
00341     }
00342 
00343     const DetectorElement* de = this->SensDetElem(*hist);
00344     if (!de) return false;
00345 
00346     // wangzhe QE calculation starts here.
00347     int pmtid = this->SensDetId(*de);
00348     DayaBay::Detector detector(pmtid);
00349     G4Track* track = step->GetTrack();
00350     double weight = track->GetWeight();
00351 
00352     // Make sure the energy deposit is caused by an optical photon.
00353     // In the future we may decide that energy deposits in the photocathode by other particles
00354     // are capable of generating a PMT signal. In that case, we will need a model for the
00355     // effects. djaffe
00356 
00357     if ( track->GetDefinition() != G4OpticalPhoton::Definition() ) 
00358       {
00359         debug() << "Track definition " << track->GetDefinition()->GetParticleName() 
00360                 << "This track is not an optical photon. abort further processing " << endreq;
00361         return false;
00362       }
00363 
00364     double qe = this->SensDetQE(preStepPoint->GetPhysicalVolume()->GetLogicalVolume(),energyDep);
00365 
00366         // Scale total efficiency upward to allow for PMT-to-PMT
00367         // efficiency variation applied in the electronics simulation.
00368         // The QEScale used here should be the inverse of the mean
00369         // PMT efficiency applied in the electronics simulation.
00370         qe *= m_qeScale;
00371       
00372         // For now, hard code "extra" decrease in efficiency for water shield PMTs to match G4dyb.
00373         if (detector.detectorId() == DetectorId::kIWS || detector.detectorId() == DetectorId::kOWS) {
00374           qe *= 0.6;
00375         }
00376     if (qe < 0.) return false;
00377     if ( m_applyPreQE ) qe/=m_preQE; 
00378     if (qe > 1.) {
00379             warning() << "Non physical QE > 1." << endreq;
00380         }
00381 
00382         // apply QE
00383         double looser = m_uni();
00384         // debug()<<"qe= "<<qe<<"   looser= "<<looser<<endreq;
00385         if (looser > qe) {
00386           //verbose() << "Failed quantum efficiency " << looser << " > " << qe << endreq;
00387           return true;
00388         }
00389     // wz QE calculation and sampling end here. PE is going to be created next.
00390 
00391     double wavelength = CLHEP::twopi*CLHEP::hbarc/energyDep;
00392 
00393 #if 1
00394     if (!pmtid) {
00395         warning () << "Dumping TouchableHistory:\n";
00396         for (int ind=0; ind<hist->GetHistoryDepth(); ++ind) {
00397             warning () << "\t (" << ind << ") " << hist->GetVolume(ind)->GetName() 
00398                     << " #" << hist->GetReplicaNumber(ind)
00399                     << " physvol @ " << (void*)hist->GetVolume(ind)
00400                     << "\n";
00401         }
00402         warning() << endreq;
00403     }
00404 
00405     verbose() << "Dumping TouchableHistory:\n";
00406     for (int ind=0; ind<hist->GetHistoryDepth(); ++ind) {
00407         verbose() << "\t (" << ind << ") " << hist->GetVolume(ind)->GetName() 
00408                   << " #" << hist->GetReplicaNumber(ind)
00409                   << " physvol @ " << (void*)hist->GetVolume(ind)
00410                   << "\n";
00411     }
00412     verbose() << endreq;
00413 #endif
00414 
00415 
00416     DayaBay::SimPmtHit* sphit = new DayaBay::SimPmtHit();
00417 
00418     // base hit
00419 
00420     // Time since event created
00421     sphit->setHitTime(preStepPoint->GetGlobalTime());
00422 
00423     //#include "G4NavigationHistory.hh"
00424 
00425     const G4AffineTransform& trans = hist->GetHistory()->GetTopTransform();
00426     const G4ThreeVector& global_pos = preStepPoint->GetPosition();
00427     G4ThreeVector pos = trans.TransformPoint(global_pos);
00428     sphit->setLocalPos(pos);
00429     sphit->setSensDetId(pmtid);
00430     
00431     // pmt hit
00432     // sphit->setDir(...);       // for now
00433     G4ThreeVector pol = trans.TransformAxis(track->GetPolarization());
00434     pol = pol.unit();
00435     G4ThreeVector dir = trans.TransformAxis(track->GetMomentum());
00436     dir = dir.unit();
00437     sphit->setPol(pol);
00438     sphit->setDir(dir);
00439     sphit->setWavelength(wavelength);
00440     sphit->setType(0);
00441     // G4cerr<<"PMT: set hit weight "<<weight<<G4endl; //gonchar
00442     sphit->setWeight(weight);
00443 
00444     // Dump info on the hit and touchable history
00445 #if 0
00446     debug() << "hit: id=" << (void*)pmtid << ", "
00447             << wavelength/CLHEP::nm << " nm, "
00448             << "pol=[" << pol[0] << "," << pol[1] << "," << pol[2] << "] "
00449             << "pos=[" << pos[0] << "," << pos[1] << "," << pos[2] << "] "
00450             << endreq;
00451 
00452     for (int ind=0; ind<hist->GetHistoryDepth(); ++ind) {
00453         verbose () << "\t (" << ind << ") " << hist->GetVolume(ind)->GetName() 
00454                 << " #" << hist->GetReplicaNumber(ind)
00455                 << " physvol @ " << (void*)hist->GetVolume(ind)
00456                 << "\n";
00457     }
00458     verbose() << endreq;
00459 #endif
00460 
00461 
00462     int trackid = track->GetTrackID();
00463     this->StoreHit(sphit,trackid);
00464 
00465     return true;
00466 }

StatusCode DsPmtSensDet::initialize (  )  [virtual]

Reimplemented from GiGaSensDetBase.

Definition at line 97 of file DsPmtSensDet.cc.

00098 {
00099     this->GiGaSensDetBase::initialize();
00100     info() << "DsPmtSensDet initialize" << endreq;
00101 
00102     m_t2de = tool<ITouchableToDetectorElement>(m_t2deName);
00103 
00104     // define collections.
00105 
00106     // first a catch-all
00107     collectionName.insert("unknown");
00108     for (int isite=0; site_ids[isite] >= 0; ++isite) {
00109         Site::Site_t site = site_ids[isite];
00110 
00111         for (int idet=0; detector_ids[idet] >= 0; ++idet) {
00112             DetectorId::DetectorId_t detid = detector_ids[idet];
00113 
00114             DayaBay::Detector det(site,detid);
00115 
00116             if (det.bogus()) continue;
00117 
00118             string name=det.detName();
00119             collectionName.insert(name.c_str());
00120             debug() << "Inserted collectionName " << collectionName.size()-1 
00121                     << " = \"" << name << "\"" 
00122                     << "(" << isite << "," << idet << "): "
00123                     << "(" << site << "," << detid << ")"
00124                     << endreq;
00125         }
00126     }
00127 
00128     // Set up a fast sim model for PMT (see ExN05).  Note, this is
00129     // kind of circular, but I find no other way to attach a fast sim
00130     // to a sens det.
00131 
00132     //const G4LogicalVolume* lv = GiGaVolumeUtils::findLVolume(m_cathodeLogVol);
00133     IGiGaGeomCnvSvc *geoSvc = 0;
00134     if (service("GiGaGeo",geoSvc,true).isFailure()) {
00135         fatal() << "Failed to get GiGa Geometry service" << endreq;
00136         return StatusCode::FAILURE;        
00137     }
00138     std::vector<std::string>::iterator volumeIter, 
00139       volumeEnd = m_cathodeLogVols.end();
00140     int PMTii = 1;
00141     char str[10] = "";
00142 
00143     for(volumeIter = m_cathodeLogVols.begin(); volumeIter != volumeEnd; 
00144         volumeIter++){
00145       std::string cathodeLogVol = *volumeIter;
00146       G4LogicalVolume* lv = geoSvc->g4LVolume( cathodeLogVol );
00147       if (!lv) {
00148         error() << "Failed to get " << cathodeLogVol << endreq;
00149         return StatusCode::FAILURE;
00150       }
00151       std::string PMT_region_name = "PMT_region";
00152       std::string PmtModel_name = "PmtModel";
00153       sprintf(str,"%d",PMTii);
00154       PMT_region_name += str;
00155       PmtModel_name += str;
00156       PMTii++;
00157 
00158       G4Region* pmt_reg = new G4Region(PMT_region_name.c_str());
00159       pmt_reg->AddRootLogicalVolume(const_cast<G4LogicalVolume*>(lv));
00160       
00161       vector<double> cuts; 
00162       cuts.push_back(1.0*mm);
00163       cuts.push_back(1.0*mm);
00164       cuts.push_back(1.0*mm);
00165       pmt_reg->SetProductionCuts(new G4ProductionCuts());
00166       pmt_reg->GetProductionCuts()->SetProductionCuts(cuts);
00167       
00168       DsPmtModel* pmtModel = new DsPmtModel(PmtModel_name.c_str(),pmt_reg);
00169       pmtModel = 0;
00170     }
00171 
00172     // set up random numbers
00173     IRndmGenSvc *rgs = 0;
00174     if (service("RndmGenSvc",rgs,true).isFailure()) {
00175         fatal() << "Failed to get random service" << endreq;
00176         return StatusCode::FAILURE;        
00177     }
00178     StatusCode sc = m_uni.initialize(rgs, Rndm::Flat(0,1));
00179     if (sc.isFailure()) {
00180         fatal() << "Failed to initialize uniform random numbers" << endreq;
00181         return StatusCode::FAILURE;
00182     }
00183 
00184     info()<<"Photons prescaling is "<<( m_applyPreQE?"on":"off" )
00185           <<". Preliminary applied efficiency is "<<m_preQE<<endreq;
00186     info()<<"Make sure the same values are set for the DsPhysConsOptical"<<endreq;
00187 
00188     return StatusCode::SUCCESS;
00189 }

StatusCode DsPmtSensDet::finalize (  )  [virtual]

Reimplemented from GiGaSensDetBase.

Definition at line 191 of file DsPmtSensDet.cc.

00192 {
00193     info() << "DsPmtSensDet finalize" << endreq;
00194     return this->GiGaSensDetBase::finalize();
00195 }

void DsPmtSensDet::StoreHit ( DayaBay::SimPmtHit hit,
int  trackid 
) [private]

Definition at line 468 of file DsPmtSensDet.cc.

00469 {
00470     int did = hit->sensDetId();
00471     DayaBay::Detector det(did);
00472     short int sdid = det.siteDetPackedData();
00473 
00474     G4DhHitCollection* hc = m_hc[sdid];
00475     if (!hc) {
00476         warning() << "Got hit with no hit collection.  ID = " << (void*)did
00477                   << " which is detector: \"" << DayaBay::Detector(did).detName()
00478                   << "\". Storing to the " << collectionName[0] << " collection"
00479                   << endreq;
00480         sdid = 0;
00481         hc = m_hc[sdid];
00482     }
00483 
00484 #if 1
00485     verbose() << "Storing hit PMT: " << (void*)did 
00486               << " from " << DayaBay::Detector(did).detName()
00487               << " in hc #"<<  sdid << " = "
00488               << hit->hitTime()/CLHEP::ns << "[ns] "
00489               << hit->localPos()/CLHEP::cm << "[cm] "
00490               << hit->wavelength()/CLHEP::nm << "[nm]"
00491               << endreq;
00492 #endif
00493 
00494     hc->insert(new G4DhHit(hit,trackid));
00495 }

const DetectorElement * DsPmtSensDet::SensDetElem ( const G4TouchableHistory &  hist  )  [private]

Definition at line 234 of file DsPmtSensDet.cc.

00235 {
00236     const IDetectorElement* idetelem = 0;
00237     int steps=0;
00238 
00239     if (!hist.GetHistoryDepth()) {
00240         error() << "DsPmtSensDet::SensDetElem given empty touchable history" << endreq;
00241         return 0;
00242     }
00243 
00244     StatusCode sc = 
00245         m_t2de->GetBestDetectorElement(&hist,m_sensorStructures,idetelem,steps);
00246     if (sc.isFailure()) {      // verbose warning
00247         warning() << "Failed to find detector element in:\n";
00248         for (size_t ind=0; ind<m_sensorStructures.size(); ++ind) {
00249             warning() << "\t\t" << m_sensorStructures[ind] << "\n";
00250         }
00251         warning() << "\tfor touchable history:\n";
00252         for (int ind=0; ind < hist.GetHistoryDepth(); ++ind) {
00253             warning() << "\t (" << ind << ") " 
00254                       << hist.GetVolume(ind)->GetName() << "\n";
00255         }
00256         warning() << endreq;
00257         return 0;
00258     }
00259     
00260     return dynamic_cast<const DetectorElement*>(idetelem);
00261 }

int DsPmtSensDet::SensDetId ( const DetectorElement de  )  [private]

Definition at line 263 of file DsPmtSensDet.cc.

00264 {
00265     const DetectorElement* detelem = &de;
00266 
00267     while (detelem) {
00268         if (detelem->params()->exists(m_idParameter)) {
00269             break;
00270         }
00271         detelem = dynamic_cast<const DetectorElement*>(detelem->parentIDetectorElement());
00272     }
00273     if (!detelem) {
00274         warning() << "Could not get PMT detector element starting from " << de << endreq;
00275         return 0;
00276     }
00277 
00278     return detelem->params()->param<int>(m_idParameter);
00279 }

double DsPmtSensDet::SensDetQE ( G4LogicalVolume *  logvol,
double  energy 
) [private]

Definition at line 280 of file DsPmtSensDet.cc.

00281 {
00282     G4Material* mat = logvol->GetMaterial();
00283     if (!mat) {
00284         warning () << "No material for " << logvol->GetName() << endreq;
00285         return -1;
00286     }
00287     
00288 
00289      G4MaterialPropertiesTable* mattab = mat->GetMaterialPropertiesTable();
00290     if (mattab) {
00291         G4MaterialPropertyVector* qevec = mattab->GetProperty(m_qeffParamName.c_str());
00292         if (qevec) {
00293         
00294           debug() << m_qeffParamName << ":("
00295                   << qevec->GetMinPhotonEnergy()/CLHEP::eV << " eV," 
00296                   << qevec->GetMinProperty() << ")-->("
00297                   << qevec->GetMaxPhotonEnergy()/CLHEP::eV << " eV," 
00298                   << qevec->GetMaxProperty() << ")"
00299                   << " particle energy is " << energy/CLHEP::eV
00300                   << endreq;
00301 
00302           return qevec->GetProperty(energy);
00303 
00304         }
00305     }
00306     else {
00307         debug () << "No material properties in " << logvol->GetName() << endreq;
00308     }
00309 
00310     int ndaught = logvol->GetNoDaughters();
00311     for (int ind=0; ind < ndaught; ++ind) {
00312         G4VPhysicalVolume* physvol = logvol->GetDaughter(ind);
00313         double qe = this->SensDetQE(physvol->GetLogicalVolume(),energy);
00314         if (qe < 0) return qe;
00315     }
00316     warning() << "All attempts failed to find " << m_qeffParamName 
00317               << " in " << logvol->GetName() << endreq;
00318     return -1;
00319 }


Member Data Documentation

std::vector<std::string> DsPmtSensDet::m_cathodeLogVols [private]

CathodeLogicalVolumes : name of logical volumes in which this sensitive detector is operating.

Definition at line 48 of file DsPmtSensDet.h.

std::vector<std::string> DsPmtSensDet::m_sensorStructures [private]

SensorStructures : names of paths in TDS in which to search for sensor detector elements using this sensitive detector.

Definition at line 52 of file DsPmtSensDet.h.

std::string DsPmtSensDet::m_idParameter [private]

PackedIdParameterName : name of user paramater of the counted detector element which holds the packed, globally unique PMT ID.

Definition at line 57 of file DsPmtSensDet.h.

std::string DsPmtSensDet::m_t2deName [private]

TouchableToDetelem : the ITouchableToDetectorElement to use to resolve sensor ID.

Definition at line 61 of file DsPmtSensDet.h.

ITouchableToDetectorElement* DsPmtSensDet::m_t2de [private]

Definition at line 62 of file DsPmtSensDet.h.

double DsPmtSensDet::m_qeScale [private]

QEScale: Upward adjustment of DetSim efficiency to allow PMT-to-PMT efficiency variation in the electronics simulation.

The value should be the inverse of the mean PMT efficiency applied in ElecSim.

Definition at line 68 of file DsPmtSensDet.h.

bool DsPmtSensDet::m_applyPreQE [private]

preQE is maximal possible PMT QE.

It's used in order to avoid generating in scintillation and Cerenkov processes unnecessary photons, which anyway will be cut out by the QE.

Definition at line 73 of file DsPmtSensDet.h.

double DsPmtSensDet::m_preQE [private]

Definition at line 74 of file DsPmtSensDet.h.

std::string DsPmtSensDet::m_qeffParamName [private]

QEffParameterName : name of user parameter in the photo cathode volume that holds the quantum efficiency tabproperty.

Definition at line 78 of file DsPmtSensDet.h.

bool DsPmtSensDet::m_preserveWeight [private]

PreserveWeight : If true, any track weight will simply be passed on in the Hit's weight data member.

O.w. it is multiplied to the QE.

Definition at line 83 of file DsPmtSensDet.h.

LocalHitCache DsPmtSensDet::m_hc [private]

Definition at line 96 of file DsPmtSensDet.h.

Rndm::Numbers DsPmtSensDet::m_uni [private]

Definition at line 98 of file DsPmtSensDet.h.


The documentation for this class was generated from the following files:
| Classes | Job Modules | Data Objects | Services | Algorithms | Tools | Packages | Directories | Tracs |

Generated on Mon Apr 11 20:53:28 2011 for DetSim by doxygen 1.4.7