Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-07-12 07:54:44

0001 //==========================================================================
0002 //  AIDA Detector description implementation 
0003 //--------------------------------------------------------------------------
0004 // Copyright (C) Organisation europeenne pour la Recherche nucleaire (CERN)
0005 // All rights reserved.
0006 //
0007 // For the licensing terms see $DD4hepINSTALL/LICENSE.
0008 // For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
0009 //
0010 // Author     : M.Frank
0011 //
0012 //==========================================================================
0013 
0014 // Framework include files
0015 #include <DD4hep/Printout.h>
0016 #include <DD4hep/Primitives.h>
0017 #include <DD4hep/InstanceCount.h>
0018 
0019 #include <DDG4/Geant4Kernel.h>
0020 #include <DDG4/Geant4Mapping.h>
0021 #include <DDG4/Geant4StepHandler.h>
0022 #include <DDG4/Geant4SensDetAction.h>
0023 #include <DDG4/Geant4VolumeManager.h>
0024 #include <DDG4/Geant4MonteCarloTruth.h>
0025 
0026 // Geant4 include files
0027 #include <G4Step.hh>
0028 #include <G4SDManager.hh>
0029 #include <G4VSensitiveDetector.hh>
0030 
0031 // C/C++ include files
0032 #include <stdexcept>
0033 
0034 #include "G4OpticalParameters.hh"
0035 #include "G4OpticalPhoton.hh"
0036 
0037 #ifdef DD4HEP_USE_GEANT4_UNITS
0038 #define MM_2_CM 1.0
0039 #else
0040 #define MM_2_CM 0.1
0041 #endif
0042 
0043 using namespace dd4hep::sim;
0044 
0045 #if 0
0046 namespace {
0047   Geant4ActionSD* _getSensitiveDetector(const std::string& name) {
0048     G4SDManager* mgr = G4SDManager::GetSDMpointer();
0049     G4VSensitiveDetector* sd = mgr->FindSensitiveDetector(name);
0050     if (0 == sd) {
0051       dd4hep::except("Geant4Sensitive", "DDG4: You requested to configure actions "
0052                      "for the sensitive detector %s,\nDDG4: which is not known to Geant4. "
0053                      "Are you sure you already converted the geometry?", name.c_str());
0054     }
0055     Geant4ActionSD* action_sd = dynamic_cast<Geant4ActionSD*>(sd);
0056     if (0 == action_sd) {
0057       throw dd4hep::except("Geant4Sensitive", "DDG4: You may only configure actions "
0058                            "for sensitive detectors of type Geant4ActionSD.\n"
0059                            "DDG4: The sensitive detector of %s is of type %s, which is incompatible.", name.c_str(),
0060                            typeName(typeid(*sd)).c_str());
0061     }
0062     return action_sd;
0063   }
0064 }
0065 #endif
0066 
0067 /// Standard action constructor
0068 Geant4ActionSD::Geant4ActionSD(const std::string& nam)
0069   : Geant4Action(0, nam) {
0070   InstanceCount::increment(this);
0071 }
0072 
0073 /// Default destructor
0074 Geant4ActionSD::~Geant4ActionSD() {
0075   InstanceCount::decrement(this);
0076 }
0077 
0078 /// Standard constructor
0079 Geant4Filter::Geant4Filter(Geant4Context* ctxt, const std::string& nam)
0080   : Geant4Action(ctxt, nam) {
0081   InstanceCount::increment(this);
0082 }
0083 
0084 /// Standard destructor
0085 Geant4Filter::~Geant4Filter() {
0086   InstanceCount::decrement(this);
0087 }
0088 
0089 /// Filter action. Return true if hits should be processed
0090 bool Geant4Filter::operator()(const G4Step*) const {
0091   return true;
0092 }
0093 
0094 /// Filter action. Return true if hits should be processed
0095 bool Geant4Filter::operator()(const Geant4FastSimSpot*) const {
0096   except("The filter action %s does not support the GFLASH/FastSim interface for Geant4.", c_name());
0097   return false;
0098 }
0099 
0100 /// Constructor. The detector element is identified by the name
0101 Geant4Sensitive::Geant4Sensitive(Geant4Context* ctxt, const std::string& nam, DetElement det, Detector& det_ref)
0102   : Geant4Action(ctxt, nam), m_detDesc(det_ref), m_detector(det)
0103 {
0104   InstanceCount::increment(this);
0105   if (!det.isValid()) {
0106     except("DDG4: Detector elemnt for %s is invalid.", nam.c_str());
0107   }
0108   declareProperty("HitCreationMode", m_hitCreationMode = SIMPLE_MODE);
0109   m_sequence     = context()->kernel().sensitiveAction(m_detector.name());
0110   m_sensitive    = m_detDesc.sensitiveDetector(det.name());
0111   m_readout      = m_sensitive.readout();
0112   m_segmentation = m_readout.segmentation();
0113 }
0114 
0115 /// Standard destructor
0116 Geant4Sensitive::~Geant4Sensitive() {
0117   m_filters(&Geant4Filter::release);
0118   m_filters.clear();
0119   InstanceCount::decrement(this);
0120 }
0121 
0122 /// Add an actor responding to all callbacks. Sequence takes ownership.
0123 void Geant4Sensitive::adoptFilter(Geant4Action* action)   {
0124   Geant4Filter* filter = dynamic_cast<Geant4Filter*>(action);
0125   adopt(filter);
0126 }
0127 
0128 /// Add an actor responding to all callbacks. Sequence takes ownership.
0129 void Geant4Sensitive::adopt(Geant4Filter* filter) {
0130   if (filter) {
0131     filter->addRef();
0132     m_filters.add(filter);
0133     return;
0134   }
0135   except("Attempt to add invalid sensitive filter!");
0136 }
0137 
0138 /// Add an actor responding to all callbacks. Sequence takes ownership.
0139 void Geant4Sensitive::adoptFilter_front(Geant4Action* action)   {
0140   Geant4Filter* filter = dynamic_cast<Geant4Filter*>(action);
0141   adopt_front(filter);
0142 }
0143 
0144 /// Add an actor responding to all callbacks. Sequence takes ownership.
0145 void Geant4Sensitive::adopt_front(Geant4Filter* filter) {
0146   if (filter) {
0147     filter->addRef();
0148     m_filters.add_front(filter);
0149     return;
0150   }
0151   except("Attempt to add invalid sensitive filter!");
0152 }
0153 
0154 /// Callback before hit processing starts. Invoke all filters.
0155 bool Geant4Sensitive::accept(const G4Step* step) const {
0156   bool (Geant4Filter::*filter)(const G4Step*) const = &Geant4Filter::operator();
0157   bool result = m_filters.filter(filter, step);
0158   return result;
0159 }
0160 
0161 /// GFLASH/FastSim interface: Callback before hit processing starts. Invoke all filters.
0162 bool Geant4Sensitive::accept(const Geant4FastSimSpot* spot) const {
0163   bool (Geant4Filter::*filter)(const Geant4FastSimSpot*) const = &Geant4Filter::operator();
0164   bool result = m_filters.filter(filter, spot);
0165   return result;
0166 }
0167 
0168 /// Access to the sensitive detector object
0169 void Geant4Sensitive::setDetector(Geant4ActionSD* sens_det) {
0170   m_sensitiveDetector = sens_det;
0171 }
0172 
0173 /// Access to the sensitive detector object
0174 Geant4ActionSD& Geant4Sensitive::detector() const {
0175   if (m_sensitiveDetector)
0176     return *m_sensitiveDetector;
0177   //m_sensitiveDetector = _getSensitiveDetector(m_detector.name());
0178   //if (  m_sensitiveDetector ) return *m_sensitiveDetector;
0179   except("DDG4: The sensitive detector for action %s was not properly configured.", name().c_str());
0180   throw std::runtime_error("Geant4Sensitive::detector");
0181 }
0182 
0183 /// Access to the hosting sequence
0184 Geant4SensDetActionSequence& Geant4Sensitive::sequence() const {
0185   return *m_sequence;
0186 }
0187 
0188 /// Access the detector desciption object
0189 dd4hep::Detector& Geant4Sensitive::detectorDescription() const {
0190   return m_detDesc;
0191 }
0192 
0193 /// Access HitCollection container names
0194 const std::string& Geant4Sensitive::hitCollectionName(std::size_t which) const {
0195   return sequence().hitCollectionName(which);
0196 }
0197 
0198 /// Retrieve the hits collection associated with this detector by its serial number
0199 Geant4HitCollection* Geant4Sensitive::collection(std::size_t which) {
0200   return sequence().collection(which);
0201 }
0202 
0203 /// Retrieve the hits collection associated with this detector by its collection identifier
0204 Geant4HitCollection* Geant4Sensitive::collectionByID(std::size_t id) {
0205   return sequence().collectionByID(id);
0206 }
0207 
0208 /// Define collections created by this sensitivie action object
0209 void Geant4Sensitive::defineCollections() {
0210 }
0211 
0212 /// Method invoked at the begining of each event.
0213 void Geant4Sensitive::begin(G4HCofThisEvent* /* HCE */) {
0214 }
0215 
0216 /// Method invoked at the end of each event.
0217 void Geant4Sensitive::end(G4HCofThisEvent* /* HCE */) {
0218 }
0219 
0220 /// G4VSensitiveDetector interface: Method for generating hit(s) using the information of G4Step object.
0221 bool Geant4Sensitive::process(const G4Step* /* step */, G4TouchableHistory* /* history */) {
0222   return false;
0223 }
0224 
0225 /// GFLASH/FastSim interface: Method for generating hit(s) using the information of the Geant4FastSimSpot object.
0226 bool Geant4Sensitive::processFastSim(const Geant4FastSimSpot* /* spot */, G4TouchableHistory* /* history */) {
0227   except("The sensitive action %s does not support the GFLASH/FastSim interface for Geant4.", c_name());
0228   return false;
0229 }
0230 
0231 /// Method is invoked if the event abortion is occured.
0232 void Geant4Sensitive::clear(G4HCofThisEvent* /* HCE */) {
0233 }
0234 
0235 /// Mark the track to be kept for MC truth propagation during hit processing
0236 void Geant4Sensitive::mark(const G4Track* track) const  {
0237   Geant4MonteCarloTruth* truth = context()->event().extension<Geant4MonteCarloTruth>(false);
0238   if ( truth ) truth->mark(track);
0239 }
0240 
0241 /// Mark the track of this step to be kept for MC truth propagation during hit processing
0242 void Geant4Sensitive::mark(const G4Step* step) const  {
0243   Geant4MonteCarloTruth* truth = context()->event().extension<Geant4MonteCarloTruth>(false);
0244   if ( truth ) truth->mark(step);
0245 }
0246 
0247 /// Returns the volumeID of the sensitive volume corresponding to the step
0248 long long int Geant4Sensitive::volumeID(const G4Step* step) {
0249   Geant4StepHandler stepH(step);
0250   Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();
0251   VolumeID id = volMgr.volumeID(stepH.preTouchable());
0252   return id;
0253 }
0254 
0255 /// Returns the volumeID of the sensitive volume corresponding to the touchable history
0256 long long int Geant4Sensitive::volumeID(const G4VTouchable* touchable) {
0257   Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();
0258   VolumeID id = volMgr.volumeID(touchable);
0259   return id;
0260 }
0261 
0262 /// Returns the cellID(volumeID+local coordinate encoding) of the sensitive volume corresponding to the step
0263 long long int Geant4Sensitive::cellID(const G4Step* step) {
0264   Geant4StepHandler h(step);
0265   Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();
0266   bool UsePostStepOnly = G4OpticalParameters::Instance() && G4OpticalParameters::Instance()->GetBoundaryInvokeSD() && (step->GetTrack()->GetDefinition() == G4OpticalPhoton::Definition());
0267   VolumeID volID = volMgr.volumeID(UsePostStepOnly? h.postTouchable() : h.preTouchable());
0268   if ( m_segmentation.isValid() )  {
0269     std::exception_ptr eptr;
0270     G4ThreeVector global = UsePostStepOnly? h.postPosG4() : 0.5 * (h.prePosG4()+h.postPosG4());
0271     G4ThreeVector local  = UsePostStepOnly? h.postTouchable()->GetHistory()->GetTopTransform().TransformPoint(global) :
0272                                             h.preTouchable()->GetHistory()->GetTopTransform().TransformPoint(global);
0273     Position loc(local.x()*MM_2_CM, local.y()*MM_2_CM, local.z()*MM_2_CM);
0274     Position glob(global.x()*MM_2_CM, global.y()*MM_2_CM, global.z()*MM_2_CM);
0275     try  {
0276       VolumeID cID = m_segmentation.cellID(loc, glob, volID);
0277       return cID;
0278     }
0279     catch(const std::exception& e)   {
0280       eptr = std::current_exception();
0281       error("cellID: failed to access segmentation for VolumeID: %016lX [%ld]  [%s]", volID, volID, e.what());
0282       error("....... G4-local:   (%f, %f, %f) G4-global:   (%f, %f, %f)",
0283         local.x(), local.y(), local.z(), global.x(), global.y(), global.z());
0284       error("....... TGeo-local: (%f, %f, %f) TGeo-global: (%f, %f, %f)",
0285         loc.x(), loc.y(), loc.z(), glob.x(), glob.y(), glob.z());
0286       error("....... Pre-step:  %s  SD: %s", h.volName(h.pre), h.sdName(h.pre).c_str());
0287       if ( h.post )
0288         error("....... Post-step: %s  SD: %s", h.volName(h.post), h.sdName(h.post).c_str());
0289       std::rethrow_exception(std::move(eptr));
0290     }
0291   }
0292   return volID;
0293 }
0294 
0295 /// Returns the cellID(volumeID+local coordinate encoding) of the sensitive volume corresponding to the touchable history
0296 long long int Geant4Sensitive::cellID(const G4VTouchable* touchable, const G4ThreeVector& global) {
0297   Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();
0298   VolumeID volID = volMgr.volumeID(touchable);
0299   if ( m_segmentation.isValid() )  {
0300     std::exception_ptr eptr;
0301     G4ThreeVector local  = touchable->GetHistory()->GetTopTransform().TransformPoint(global);
0302     Position loc (local.x()*MM_2_CM, local.y()*MM_2_CM, local.z()*MM_2_CM);
0303     Position glob(global.x()*MM_2_CM, global.y()*MM_2_CM, global.z()*MM_2_CM);
0304     try  {
0305       VolumeID cID = m_segmentation.cellID(loc, glob, volID);
0306       return cID;
0307     }
0308     catch(const std::exception& e)   {
0309       auto* pvol = touchable->GetVolume();
0310       auto* vol = pvol->GetLogicalVolume();
0311       auto* sd = vol->GetSensitiveDetector();
0312       eptr = std::current_exception();
0313       error("cellID: failed to access segmentation for VolumeID: %016lX [%ld]  [%s]", volID, volID, e.what());
0314       error("....... G4-local:   (%f, %f, %f) G4-global:   (%f, %f, %f)",
0315         local.x(), local.y(), local.z(), global.x(), global.y(), global.z());
0316       error("....... TGeo-local: (%f, %f, %f) TGeo-global: (%f, %f, %f)",
0317         loc.x(), loc.y(), loc.z(), glob.x(), glob.y(), glob.z());
0318       error("....... Touchable:  %s  SD: %s", vol->GetName().c_str(), sd ? sd->GetName().c_str() : "???");
0319       std::rethrow_exception(std::move(eptr));
0320     }
0321   }
0322   return volID;
0323 }
0324 
0325 /// Standard constructor
0326 Geant4SensDetActionSequence::Geant4SensDetActionSequence(Geant4Context* ctxt, const std::string& nam)
0327   : Geant4Action(ctxt, nam), m_hce(0), m_detector(0)
0328 {
0329   m_needsControl = true;
0330   context()->sensitiveActions().insert(name(), this);
0331   /// Update the sensitive detector type, so that the proper instance is created
0332   m_sensitive = context()->detectorDescription().sensitiveDetector(nam);
0333   m_sensitiveType = m_sensitive.type();
0334   InstanceCount::increment(this);
0335 }
0336 
0337 /// Default destructor
0338 Geant4SensDetActionSequence::~Geant4SensDetActionSequence() {
0339   m_filters(&Geant4Filter::release);
0340   m_actors(&Geant4Sensitive::release);
0341   m_filters.clear();
0342   m_actors.clear();
0343   InstanceCount::decrement(this);
0344 }
0345 
0346 /// Set or update client context
0347 void Geant4SensDetActionSequence::updateContext(Geant4Context* ctxt)    {
0348   m_context = ctxt;
0349   m_actors.updateContext(ctxt);
0350   m_filters.updateContext(ctxt);
0351 }
0352 
0353 /// Add an actor responding to all callbacks. Sequence takes ownership.
0354 void Geant4SensDetActionSequence::adoptFilter(Geant4Action* action)   {
0355   Geant4Filter* filter = dynamic_cast<Geant4Filter*>(action);
0356   adopt(filter);
0357 }
0358 
0359 /// Add an actor responding to all callbacks
0360 void Geant4SensDetActionSequence::adopt(Geant4Sensitive* sensitive) {
0361   if (sensitive) {
0362     sensitive->addRef();
0363     m_actors.add(sensitive);
0364     return;
0365   }
0366   except("Attempt to add invalid sensitive actor!");
0367 }
0368 
0369 /// Add an actor responding to all callbacks. Sequence takes ownership.
0370 void Geant4SensDetActionSequence::adopt(Geant4Filter* filter) {
0371   if (filter) {
0372     filter->addRef();
0373     m_filters.add(filter);
0374     return;
0375   }
0376   except("Attempt to add invalid sensitive filter!");
0377 }
0378 
0379 /// Initialize the usage of a hit collection. Returns the collection identifier
0380 std::size_t Geant4SensDetActionSequence::defineCollection(Geant4Sensitive* owner, const std::string& collection_name, create_t func) {
0381   m_collections.emplace_back(collection_name, make_pair(owner,func));
0382   return m_collections.size() - 1;
0383 }
0384 
0385 /// Called at construction time of the sensitive detector to declare all hit collections
0386 std::size_t Geant4SensDetActionSequence::Geant4SensDetActionSequence::defineCollections(Geant4ActionSD* sens_det) {
0387   std::size_t count = 0;
0388   m_detector = sens_det;
0389   m_actors(&Geant4Sensitive::setDetector, sens_det);
0390   m_actors(&Geant4Sensitive::defineCollections);
0391   for (HitCollections::const_iterator i = m_collections.begin(); i != m_collections.end(); ++i) {
0392     sens_det->defineCollection((*i).first);
0393     ++count;
0394   }
0395   return count;
0396 }
0397 
0398 /// Access HitCollection container names
0399 const std::string& Geant4SensDetActionSequence::hitCollectionName(std::size_t which) const {
0400   if (which < m_collections.size()) {
0401     return m_collections[which].first;
0402   }
0403   static std::string blank = "";
0404   except("The collection name index for subdetector %s is out of range!", c_name());
0405   return blank;
0406 }
0407 
0408 /// Retrieve the hits collection associated with this detector by its serial number
0409 Geant4HitCollection* Geant4SensDetActionSequence::collection(std::size_t which) const {
0410   if (which < m_collections.size()) {
0411     int hc_id = m_detector->GetCollectionID(which);
0412     Geant4HitCollection* c = (Geant4HitCollection*) m_hce->GetHC(hc_id);
0413     if (c)
0414       return c;
0415     except("The collection index for subdetector %s is wrong!", c_name());
0416   }
0417   except("The collection name index for subdetector %s is out of range!", c_name());
0418   return 0;
0419 }
0420 
0421 /// Retrieve the hits collection associated with this detector by its collection identifier
0422 Geant4HitCollection* Geant4SensDetActionSequence::collectionByID(std::size_t id) const {
0423   Geant4HitCollection* c = (Geant4HitCollection*) m_hce->GetHC(id);
0424   if (c)
0425     return c;
0426   except("The collection index for subdetector %s is wrong!", c_name());
0427   return 0;
0428 }
0429 
0430 /// Callback before hit processing starts. Invoke all filters.
0431 bool Geant4SensDetActionSequence::accept(const G4Step* step) const {
0432   bool (Geant4Filter::*filter)(const G4Step*) const = &Geant4Filter::operator();
0433   bool result = m_filters.filter(filter, step);
0434   return result;
0435 }
0436 
0437 /// Callback before hit processing starts. Invoke all filters.
0438 bool Geant4SensDetActionSequence::accept(const Geant4FastSimSpot* spot) const {
0439   bool (Geant4Filter::*filter)(const Geant4FastSimSpot*) const = &Geant4Filter::operator();
0440   bool result = m_filters.filter(filter, spot);
0441   return result;
0442 }
0443 
0444 /// G4VSensitiveDetector interface: Method for generating hit(s) using the information of G4Step object.
0445 bool Geant4SensDetActionSequence::process(const G4Step* step, G4TouchableHistory* history) {
0446   bool result = false;
0447   for (Geant4Sensitive* sensitive : m_actors)  {
0448     if ( sensitive->accept(step) )
0449       result |= sensitive->process(step, history);
0450   }
0451   m_process(step, history);
0452   return result;
0453 }
0454 
0455 /// GFLASH/FastSim interface: Method for generating hit(s) using the information of the Geant4FastSimSpot object.
0456 bool Geant4SensDetActionSequence::processFastSim(const Geant4FastSimSpot* spot, G4TouchableHistory* history)  {
0457   bool result = false;
0458   for (Geant4Sensitive* sensitive : m_actors)  {
0459     if ( sensitive->accept(spot) )
0460       result |= sensitive->processFastSim(spot, history);
0461   }
0462   m_process(spot, history);
0463   return result;
0464 }
0465 
0466 /** G4VSensitiveDetector interface: Method invoked at the begining of each event.
0467  *  The hits collection(s) created by this sensitive detector must
0468  *  be set to the G4HCofThisEvent object at one of these two methods.
0469  */
0470 void Geant4SensDetActionSequence::begin(G4HCofThisEvent* hce) {
0471   m_hce = hce;
0472   for (std::size_t count = 0; count < m_collections.size(); ++count) {
0473     const HitCollection& cr = m_collections[count];
0474     Geant4HitCollection* col = (*cr.second.second)(name(), cr.first, cr.second.first);
0475     int id = m_detector->GetCollectionID(count);
0476     m_hce->AddHitsCollection(id, col);
0477   }
0478   m_actors(&Geant4Sensitive::begin, m_hce);
0479   m_begin (m_hce);
0480 }
0481 
0482 /// G4VSensitiveDetector interface: Method invoked at the end of each event.
0483 void Geant4SensDetActionSequence::end(G4HCofThisEvent* hce) {
0484   m_end(hce);
0485   m_actors(&Geant4Sensitive::end, hce);
0486   // G4HCofThisEvent must be availible until end-event. m_hce = 0;
0487 }
0488 
0489 /// G4VSensitiveDetector interface: Method invoked if the event was aborted.
0490 /** Hits collections created but not being set to G4HCofThisEvent
0491  *  at the event should be deleted.
0492  *  Collection(s) which have already set to G4HCofThisEvent
0493  *  will be deleted automatically.
0494  */
0495 void Geant4SensDetActionSequence::clear() {
0496   m_clear (m_hce);
0497   m_actors(&Geant4Sensitive::clear, m_hce);
0498 }
0499 
0500 /// Default destructor
0501 Geant4SensDetSequences::~Geant4SensDetSequences() {
0502   detail::releaseObjects(m_sequences);
0503   m_sequences.clear();
0504 }
0505 
0506 /// Access sequence member by name
0507 Geant4SensDetActionSequence* Geant4SensDetSequences::operator[](const std::string& nam) const {
0508   std::string n = "SD_Seq_" + nam;
0509   Members::const_iterator i = m_sequences.find(n);
0510   if (i != m_sequences.end())
0511     return (*i).second;
0512   except("Attempt to access undefined SensDetActionSequence: %s ", nam.c_str());
0513   return nullptr;
0514 }
0515 
0516 /// Access sequence member by name
0517 Geant4SensDetActionSequence* Geant4SensDetSequences::find(const std::string& name) const {
0518   std::string nam = "SD_Seq_" + name;
0519   Members::const_iterator i = m_sequences.find(nam);
0520   if (i != m_sequences.end())
0521     return (*i).second;
0522   return 0;
0523 }
0524 
0525 /// Insert sequence member
0526 void Geant4SensDetSequences::insert(const std::string& name, Geant4SensDetActionSequence* seq) {
0527   if (seq) {
0528     std::string nam = "SD_Seq_" + name;
0529     seq->addRef();
0530     m_sequences[nam] = seq;
0531     return;
0532   }
0533   except("Attempt to add invalid sensitive sequence with name:%s", name.c_str());
0534 }
0535 
0536 /// Clear the sequence list
0537 void Geant4SensDetSequences::clear()   {
0538   m_sequences.clear();
0539 }