Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-03-13 08:20:08

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 <DDG4/Geant4SensDetAction.inl>
0016 #include <DDG4/Geant4FastSimHandler.h>
0017 #include <DDG4/Geant4EventAction.h>
0018 #include <G4OpticalPhoton.hh>
0019 #include <G4VProcess.hh>
0020 
0021 
0022 /// Namespace for the AIDA detector description toolkit
0023 namespace dd4hep {
0024   
0025   /// Namespace for the Geant4 based simulation part of the AIDA detector description toolkit
0026   namespace sim   {
0027 
0028     namespace {
0029       struct Geant4VoidSensitive {};
0030 
0031       /// Common code to handle the creation of a calorimeter hit.
0032       template <class HANDLER>
0033       void handleCalorimeterHit (VolumeID cell,
0034                                  const HitContribution& contrib,
0035                                  Geant4HitCollection& coll,
0036                                  const HANDLER& h,
0037                                  const Geant4Sensitive& sd,
0038                                  const Segmentation& segmentation)
0039       {
0040         typedef Geant4Calorimeter::Hit Hit;
0041         Hit* hit = coll.findByKey<Hit>(cell);
0042         if ( !hit ) {
0043           DDSegmentation::Vector3D pos = segmentation.position(cell);
0044           Position global;
0045 
0046           // Convert the position relative to the local readout volume
0047           // to a global position.
0048           if (!segmentation.cellsSpanVolumes()) {
0049             global = h.localToGlobal(pos);
0050           }
0051           else {
0052             // The segmentation can gang together multiple volumes.
0053             // In this case, we can't use the transformation we get from
0054             // the step --- the volume that actually contains the hit
0055             // may not be the same volume that the segmentation uses
0056             // for the local coordinate system.  We need to get the
0057             // actual volID used from the segmentation and then look
0058             // it up the volume manager to get the proper transformation.
0059             VolumeID volID = segmentation.volumeID(cell);
0060             VolumeManager vman = VolumeManager::getVolumeManager(sd.detectorDescription());
0061             VolumeManagerContext* vc = vman.lookupContext(volID);
0062             // explicit unit conversion; h.localToGlobal does it internally already
0063             global = vc->localToWorld(Position(pos)) / dd4hep::mm;
0064           }
0065           hit = new Hit(global);
0066           hit->cellID = cell;
0067           coll.add(cell, hit);
0068           Geant4TouchableHandler handler(h.touchable());
0069           sd.printM2("%s> CREATE hit with deposit:%e MeV  Pos:%8.2f %8.2f %8.2f  %s  [%s]",
0070                      sd.c_name(),contrib.deposit,pos.X,pos.Y,pos.Z,handler.path().c_str(),
0071                      coll.GetName().c_str());
0072           if ( 0 == hit->cellID )  { // for debugging only!
0073             hit->cellID = cell;
0074             sd.except("+++ Invalid CELL ID for hit!");
0075           }
0076         }
0077         hit->truth.emplace_back(contrib);
0078         hit->energyDeposit += contrib.deposit;
0079       }
0080     }
0081 
0082     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0083     //               Geant4SensitiveAction<Geant4VoidSensitive>
0084     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0085     /** \addtogroup Geant4SDActionPlugin
0086      *
0087      * @{
0088      * \package Geant4VoidSensitiveAction
0089      * \brief Void Sensitive detector action to skip the processing of a detector
0090      *        without changing the entire DDG4 setup.
0091      *
0092      * @}
0093      */
0094 
0095     /// Define collections created by this sensitivie action object
0096     template <> void Geant4SensitiveAction<Geant4VoidSensitive>::defineCollections()    {
0097       m_collectionID = -1;
0098     }
0099 
0100     /// G4VSensitiveDetector interface: Method for generating hit(s) using the information of G4Step object.
0101     template <> bool
0102     Geant4SensitiveAction<Geant4VoidSensitive>::process(const G4Step*       /* step */,
0103                             G4TouchableHistory* /* hist */) {
0104       return true;
0105     }
0106     
0107     /// GFLASH/FastSim interface: Method for generating hit(s) using the information of G4Step object.
0108     template <> bool
0109     Geant4SensitiveAction<Geant4VoidSensitive>::processFastSim(const Geant4FastSimSpot* /* spot */,
0110                                    G4TouchableHistory*      /* hist */) {
0111       return true;
0112     }
0113     typedef Geant4SensitiveAction<Geant4VoidSensitive> Geant4VoidSensitiveAction;
0114 
0115     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0116     //               Geant4SensitiveAction<Geant4Tracker>
0117     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0118     /** \addtogroup Geant4SDActionPlugin
0119      *
0120      * @{
0121      * \package Geant4TrackerAction
0122      * \brief Sensitive detector meant for tracking detectors, will produce one hit per step
0123      *
0124      * @}
0125      */
0126 
0127     /// Define collections created by this sensitivie action object
0128     template <> void Geant4SensitiveAction<Geant4Tracker>::defineCollections()    {
0129       m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
0130     }
0131 
0132     /// Method for generating hit(s) using the information of G4Step object.
0133     template <> bool
0134     Geant4SensitiveAction<Geant4Tracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
0135       typedef Geant4Tracker::Hit Hit;
0136       // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
0137       //       2) The energy deposit is the difference between incoming and outcoming energy.
0138       Geant4StepHandler h(step);
0139       auto contrib = Hit::extractContribution(step);
0140       Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
0141       double    hit_deposit  = h.deposit();
0142       Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
0143 
0144       hit->cellID = cellID(step);
0145       if ( 0 == hit->cellID )  {
0146         hit->cellID = volumeID(step);
0147         except("+++ Invalid CELL ID for hit!");
0148       }
0149       collection(m_collectionID)->add(hit);
0150       mark(h.track);
0151       print("Hit with deposit:%f  Pos:%f %f %f ID=%016X",
0152             hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
0153       Geant4TouchableHandler handler(step);
0154       print("    Geant4 path:%s",handler.path().c_str());
0155       return true;
0156     }
0157     
0158     /// GFlash/FastSim interface: Method for generating hit(s) using the information of Geant4FastSimSpot object.
0159     template <> bool
0160     Geant4SensitiveAction<Geant4Tracker>::processFastSim(const Geant4FastSimSpot* spot,
0161                              G4TouchableHistory* /* hist */)
0162     {
0163       typedef Geant4Tracker::Hit Hit;
0164       Geant4FastSimHandler h(spot);
0165       auto contrib = Hit::extractContribution(spot);
0166       Hit* hit = new Hit(contrib, h.momentum(), h.deposit());
0167 
0168       hit->cellID  = cellID(h.touchable(), h.avgPositionG4());
0169       if ( 0 == hit->cellID )  {
0170         hit->cellID = volumeID(h.touchable());
0171         except("+++ Invalid CELL ID for hit!");
0172       }
0173       collection(m_collectionID)->add(hit);
0174       mark(h.track);
0175       print("Hit with deposit:%f  Pos:%f %f %f ID=%016X",
0176             hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
0177       Geant4TouchableHandler handler(h.touchable());
0178       print("    Geant4 path:%s",handler.path().c_str());
0179       return true;
0180     }
0181 
0182     typedef Geant4SensitiveAction<Geant4Tracker> Geant4TrackerAction;
0183 
0184     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0185     //               Geant4SensitiveAction<Geant4PhotonCounter>
0186     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0187     /** \addtogroup Geant4SDActionPlugin
0188      *
0189      * @{
0190      * \package Geant4OpticalTrackerAction
0191      * \brief Sensitive detector meant for photon detectors, will produce one hit per step
0192      *        for regular particles, but absorb optical photons fully on their first hit
0193      *
0194      * @}
0195      */
0196 
0197     /// Helper class to define properties of optical trackers. UNTESTED
0198     struct Geant4OpticalTracker {};
0199 
0200     /// Define collections created by this sensitivie action object
0201     template <> void Geant4SensitiveAction<Geant4OpticalTracker>::defineCollections()    {
0202       m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
0203     }
0204 
0205     /// Method for generating hit(s) using the information of G4Step object.
0206     template <> bool
0207     Geant4SensitiveAction<Geant4OpticalTracker>::process(const G4Step* step,G4TouchableHistory* /* hist */) {
0208       // Note: 1) We store in the hit the hit-direction, which is not the same as the track direction.
0209       //       2) The energy deposit is the track momentum
0210       typedef Geant4Tracker::Hit Hit;
0211       Geant4StepHandler h(step);
0212       auto      contrib = Hit::extractContribution(step);
0213       Direction hit_momentum = 0.5 * (h.preMom() + h.postMom());
0214       double    hit_deposit  = contrib.deposit;
0215       Hit* hit = new Hit(contrib, hit_momentum, hit_deposit);
0216 
0217       if (h.trackDef() == G4OpticalPhoton::OpticalPhotonDefinition()) {
0218         step->GetTrack()->SetTrackStatus(fStopAndKill);
0219       }
0220       hit->cellID = cellID(step);
0221       if ( 0 == hit->cellID )  {
0222         hit->cellID      = volumeID( step ) ;
0223         except("+++ Invalid CELL ID for hit!");
0224       }
0225       collection(m_collectionID)->add(hit);
0226       mark(h.track);
0227       print("Hit with deposit:%f  Pos:%f %f %f ID=%016X",
0228             hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
0229       Geant4TouchableHandler handler(step);
0230       print("    Geant4 path:%s",handler.path().c_str());
0231       return true;
0232     }
0233 
0234     /// GFlash/FastSim interface: Method for generating hit(s) using the information of G4Step object.
0235     template <> bool
0236     Geant4SensitiveAction<Geant4OpticalTracker>::processFastSim(const Geant4FastSimSpot* spot,
0237                                     G4TouchableHistory* /* hist */)
0238     {
0239       typedef Geant4Tracker::Hit Hit;
0240       Geant4FastSimHandler h(spot);
0241       auto contrib = Hit::extractContribution(spot);
0242       Hit* hit = new Hit(contrib, h.momentum(), contrib.deposit);
0243 
0244       hit->cellID  = cellID(h.touchable(), h.avgPositionG4());
0245       if ( 0 == hit->cellID )  {
0246         hit->cellID      = volumeID(h.touchable());
0247         except("+++ Invalid CELL ID for hit!");
0248       }
0249       collection(m_collectionID)->add(hit);
0250       mark(h.track);
0251       print("Hit with deposit:%f  Pos:%f %f %f ID=%016X",
0252             hit->energyDeposit,hit->position.X(),hit->position.Y(),hit->position.Z(),(void*)hit->cellID);
0253       Geant4TouchableHandler handler(h.touchable());
0254       print("    Geant4 path:%s",handler.path().c_str());
0255       return true;
0256     }
0257     typedef Geant4SensitiveAction<Geant4OpticalTracker> Geant4OpticalTrackerAction;
0258 
0259     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0260     //               Geant4SensitiveAction<Calorimeter>
0261     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0262     /** \addtogroup Geant4SDActionPlugin
0263      *
0264      * @{
0265      * \package Geant4CalorimeterAction
0266      *
0267      * \brief Sensitive detector meant for calorimeters
0268 
0269      *
0270      * @}
0271      */
0272     /// Define collections created by this sensitivie action object
0273     template <> void Geant4SensitiveAction<Geant4Calorimeter>::defineCollections() {
0274       m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
0275     }
0276 
0277     /// Method for generating hit(s) using the information of G4Step object.
0278     template <> bool
0279     Geant4SensitiveAction<Geant4Calorimeter>::process(const G4Step* step,G4TouchableHistory*) {
0280       typedef Geant4Calorimeter::Hit Hit;
0281       Geant4StepHandler    h(step);
0282       HitContribution      contrib = Hit::extractContribution(step);
0283       Geant4HitCollection* coll    = collection(m_collectionID);
0284       VolumeID cell = 0;
0285 
0286       try {
0287         cell = cellID(step);
0288       } catch(std::runtime_error &e) {
0289         std::stringstream out;
0290         out << std::setprecision(20) << std::scientific;
0291         out << "ERROR: " << e.what()  << std::endl;
0292         out << "Position: "
0293             << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
0294             << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
0295             << std::endl;
0296         out << "Momentum: "
0297             << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum()  << ") "
0298             << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
0299             << std::endl;
0300         std::cout << out.str();
0301         return true;
0302       }
0303 
0304       handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
0305       mark(h.track);
0306       return true;
0307     }
0308     /// GFlash/FastSim interface: Method for generating hit(s) using the information of Geant4FastSimSpot object.
0309     template <> bool
0310     Geant4SensitiveAction<Geant4Calorimeter>::processFastSim(const Geant4FastSimSpot* spot,
0311                                  G4TouchableHistory* /* hist */)
0312     {
0313       typedef Geant4Calorimeter::Hit Hit;
0314       Geant4FastSimHandler h(spot);
0315       HitContribution      contrib = Hit::extractContribution(spot);
0316       Geant4HitCollection* coll    = collection(m_collectionID);
0317       VolumeID cell = 0;
0318       
0319       try {
0320         cell = cellID(h.touchable(), h.avgPositionG4());
0321       } catch(std::runtime_error &e) {
0322         std::stringstream out;
0323         out << std::setprecision(20) << std::scientific;
0324         out << "ERROR: " << e.what()  << std::endl;
0325         out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
0326         out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
0327         std::cout << out.str();
0328         return true;
0329       }
0330       handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
0331       mark(h.track);
0332       return true;
0333     }
0334 
0335     typedef Geant4SensitiveAction<Geant4Calorimeter> Geant4CalorimeterAction;
0336 
0337     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0338     //               Geant4SensitiveAction<OpticalCalorimeter>
0339     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0340 
0341     /**
0342      *  \author  M.Frank
0343      *  \version 1.0
0344      *  \addtogroup Geant4SDActionPlugin
0345      *
0346      * @{
0347      * \package Geant4OpticalCalorimeterAction
0348      *
0349      * \brief Sensitive detector meant for optical calorimeters
0350      *
0351      * @}
0352      */
0353     /// Helper class to define properties of optical calorimeters. UNTESTED
0354     struct Geant4OpticalCalorimeter {};
0355 
0356     /// Define collections created by this sensitivie action object
0357     template <> void Geant4SensitiveAction<Geant4OpticalCalorimeter>::defineCollections() {
0358       m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
0359     }
0360     /// Method for generating hit(s) using the information of G4Step object.
0361     template <> bool
0362     Geant4SensitiveAction<Geant4OpticalCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
0363       G4Track * track =  step->GetTrack();
0364       // check that particle is optical photon:
0365       if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() )  {
0366         return false;
0367       }
0368       else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov")  {
0369         track->SetTrackStatus(fStopAndKill);
0370         return false;
0371       }
0372       else {
0373         typedef Geant4Calorimeter::Hit Hit;
0374         Geant4StepHandler h(step);
0375         Geant4HitCollection*  coll    = collection(m_collectionID);
0376         HitContribution       contrib = Hit::extractContribution(step);
0377         Position              pos     = h.prePos();
0378         Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
0379         if ( !hit ) {
0380           hit = new Hit(pos);
0381           hit->cellID = volumeID(step);
0382           coll->add(hit);
0383           if ( 0 == hit->cellID )  {
0384             hit->cellID = volumeID(step);
0385             except("+++ Invalid CELL ID for hit!");
0386           }
0387         }
0388         hit->energyDeposit += contrib.deposit;
0389         hit->truth.emplace_back(contrib);
0390         track->SetTrackStatus(fStopAndKill); // don't step photon any further
0391         mark(h.track);
0392         return true;
0393       }
0394     }
0395 
0396     /// GFlash/FastSim interface: Method for generating hit(s) using the information of Geant4FastSimSpot object.
0397     template <> bool
0398     Geant4SensitiveAction<Geant4OpticalCalorimeter>::processFastSim(const Geant4FastSimSpot* spot,
0399                                     G4TouchableHistory* /* hist */)
0400     {
0401       typedef Geant4Calorimeter::Hit Hit;
0402       Geant4FastSimHandler   h(spot);
0403       const G4Track* track = h.track;
0404       if( track->GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition() )  {
0405         return false;
0406       }
0407       else if ( track->GetCreatorProcess()->G4VProcess::GetProcessName() != "Cerenkov")  {
0408         return false;
0409       }
0410       else {
0411         Geant4HitCollection* coll = collection(m_collectionID);
0412         HitContribution   contrib = Hit::extractContribution(spot);
0413         Position          pos     = h.avgPosition();
0414         Hit* hit = coll->find<Hit>(PositionCompare<Hit,Position>(pos));
0415         if ( !hit ) {
0416           hit = new Hit(pos);
0417           hit->cellID = volumeID(h.touchable());
0418           coll->add(hit);
0419           if ( 0 == hit->cellID )  {
0420             hit->cellID = volumeID(h.touchable());
0421             except("+++ Invalid CELL ID for hit!");
0422           }
0423         }
0424         hit->energyDeposit += contrib.deposit;
0425         hit->truth.emplace_back(contrib);
0426         mark(h.track);
0427         return true;
0428       }
0429     }
0430     typedef Geant4SensitiveAction<Geant4OpticalCalorimeter>  Geant4OpticalCalorimeterAction;
0431 
0432 
0433     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0434     //               Geant4SensitiveAction<ScintillatorCalorimeter>
0435     //               For scintillator with Geant4 BirksLaw effect
0436     // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
0437 
0438 
0439     /* \addtogroup Geant4SDActionPlugin
0440      *
0441      * @{
0442      * \package Geant4ScintillatorCalorimeterAction
0443      * \brief Sensitive detector meant for scintillator calorimeters
0444      *
0445      * This sensitive action will apply Birks' law to the energy deposits
0446      *
0447      * @}
0448      */
0449     /// Class to implement the standard sensitive detector for scintillator calorimeters
0450     struct Geant4ScintillatorCalorimeter {};
0451 
0452     /// Define collections created by this sensitivie action object
0453     template <> void Geant4SensitiveAction<Geant4ScintillatorCalorimeter>::defineCollections() {
0454       m_collectionID = declareReadoutFilteredCollection<Geant4Calorimeter::Hit>();
0455     }
0456     /// Method for generating hit(s) using the information of G4Step object.
0457     template <> bool
0458     Geant4SensitiveAction<Geant4ScintillatorCalorimeter>::process(const G4Step* step,G4TouchableHistory*) {
0459       typedef Geant4Calorimeter::Hit Hit;
0460       Geant4StepHandler    h(step);
0461       HitContribution      contrib = Hit::extractContribution(step,true);
0462       Geant4HitCollection* coll    = collection(m_collectionID);
0463       VolumeID cell = 0;
0464       try {
0465         cell = cellID(step);
0466       } catch(std::runtime_error &e) {
0467         std::stringstream out;
0468         out << std::setprecision(20) << std::scientific;
0469         out << "ERROR: " << e.what()  << std::endl;
0470         out << "Position: "
0471             << "Pre (" << std::setw(24) << step->GetPreStepPoint()->GetPosition() << ") "
0472             << "Post (" << std::setw(24) << step->GetPostStepPoint()->GetPosition() << ") "
0473             << std::endl;
0474         out << "Momentum: "
0475             << " Pre (" <<std::setw(24) << step->GetPreStepPoint() ->GetMomentum()  << ") "
0476             << " Post (" <<std::setw(24) << step->GetPostStepPoint()->GetMomentum() << ") "
0477             << std::endl;
0478         std::cout << out.str();
0479         return true;
0480       }
0481       handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
0482       mark(h.track);
0483       return true;
0484     }
0485 
0486     /// GFlash/FastSim interface: Method for generating hit(s) using the information of Geant4FastSimSpot object.
0487     template <> bool
0488     Geant4SensitiveAction<Geant4ScintillatorCalorimeter>::processFastSim(const Geant4FastSimSpot* spot,
0489                                      G4TouchableHistory* /* hist */)
0490     {
0491       typedef Geant4Calorimeter::Hit Hit;
0492       Geant4FastSimHandler h(spot);
0493       HitContribution         contrib = Hit::extractContribution(spot);
0494       Geant4HitCollection*    coll    = collection(m_collectionID);
0495       VolumeID cell = 0;
0496       
0497       try {
0498         cell = cellID(h.touchable(), h.avgPositionG4());
0499       } catch(std::runtime_error &e) {
0500         std::stringstream out;
0501         out << std::setprecision(20) << std::scientific;
0502         out << "ERROR: " << e.what()  << std::endl;
0503         out << "Position: (" << std::setw(24) << h.avgPositionG4() << ") " << std::endl;
0504         out << "Momentum: (" << std::setw(24) << h.momentumG4() << ") " << std::endl;
0505         std::cout << out.str();
0506         return true;
0507       }
0508       handleCalorimeterHit(cell, contrib, *coll, h, *this, m_segmentation);
0509       mark(h.track);
0510       return true;
0511     }
0512 
0513     typedef Geant4SensitiveAction<Geant4ScintillatorCalorimeter> Geant4ScintillatorCalorimeterAction;
0514 
0515     /**
0516      *
0517      * \addtogroup Geant4SDActionPlugin
0518      * @{
0519      * \package Geant4TrackerCombineAction
0520      *
0521      * \brief Sensitive detector meant for tracking detectors will combine
0522      * multiple steps of the same track in the same sensitive volume into a
0523      * single hit
0524 
0525      *
0526      * @}
0527      */
0528     /// Geant4 sensitive detector combining all deposits of one G4Track within one sensitive element.
0529     /**
0530      *  Geant4SensitiveAction<TrackerCombine>
0531      *
0532      *
0533      *  \author  M.Frank
0534      *  \version 1.0
0535      *  \ingroup DD4HEP_SIMULATION
0536      */
0537     struct TrackerCombine {
0538       Geant4Tracker::Hit pre, post;
0539       Position           mean_pos;
0540       Geant4Sensitive*   sensitive { nullptr };
0541       G4VPhysicalVolume* firstSpotVolume  { nullptr };
0542       double             mean_time {  0e0 };
0543       double             e_cut     {  0e0 };
0544       int                current   {   -1 };
0545       int                combined  {    0 };
0546       long long int      cell      {    0 };
0547 
0548       TrackerCombine() : pre(), post()   {
0549       }
0550 
0551       /// Start a new hit
0552       void start_collecting(const G4Track* track)   {
0553         pre.truth.deposit = 0.0;
0554         current = pre.truth.trackID;
0555         sensitive->mark(track);
0556         mean_pos.SetXYZ(0,0,0);
0557         mean_time = 0;
0558         post.copyFrom(pre);
0559         combined = 0;
0560         cell = 0;
0561       }
0562       void start(const G4Step* step, const G4StepPoint* point)   {
0563         pre.storePoint(step,point);
0564         start_collecting(step->GetTrack());
0565         firstSpotVolume = step->GetPreStepPoint()->GetTouchableHandle()->GetVolume();
0566       }
0567       void start(const Geant4FastSimSpot* spot)   {
0568         pre.storePoint(spot);
0569         start_collecting(spot->primary);
0570         firstSpotVolume = spot->volume();
0571       }
0572 
0573       /// Update energy and track information during hit info accumulation
0574       void update_collected_hit(const G4VTouchable* h, G4ThreeVector&& global)    {
0575         pre.truth.deposit += post.truth.deposit;
0576         mean_pos.SetX(mean_pos.x()+post.position.x()*post.truth.deposit);
0577         mean_pos.SetY(mean_pos.y()+post.position.y()*post.truth.deposit);
0578         mean_pos.SetZ(mean_pos.z()+post.position.z()*post.truth.deposit);
0579         mean_time += post.truth.time*post.truth.deposit;
0580         if ( 0 == cell )   {
0581           cell = sensitive->cellID(h, global);
0582           if ( 0 == cell )  {
0583             cell = sensitive->volumeID(h) ;
0584             sensitive->except("+++ Invalid CELL ID for hit!");
0585           }
0586         }
0587         ++combined;
0588       }
0589       void update(const Geant4StepHandler& h) {
0590         post.storePoint(h.step, h.post);
0591         update_collected_hit(h.preTouchable(), h.avgPositionG4()); // Compute cellID
0592       }
0593       void update(const Geant4FastSimHandler& h)   {
0594         post.storePoint(h.spot);
0595         update_collected_hit(h.touchable(), h.avgPositionG4());       // Compute cellID
0596       }
0597 
0598       /// Clear collected information and restart for new hit
0599       void clear()  {
0600         mean_pos.SetXYZ(0,0,0);
0601         mean_time = 0;
0602         post.clear();
0603         pre.clear();
0604         current = -1;
0605         combined = 0;
0606         cell = 0;
0607         firstSpotVolume = nullptr;
0608       }
0609 
0610       /// Helper function to decide if the hit has to be extracted and saved in the collection
0611       bool mustSaveTrack(const G4Track* tr)  const   {
0612         return current > 0 && current != tr->GetTrackID();
0613       }
0614 
0615       /// Extract hit information and add the created hit to the collection
0616       void extractHit(Geant4HitCollection* collection)   {
0617         if ( current == -1 ) {
0618           return;
0619         }
0620         double   depo = pre.truth.deposit;
0621         double   time = depo != 0 ? mean_time / depo : mean_time;
0622         Position pos  = depo != 0 ? mean_pos  / depo : mean_pos;
0623         Momentum mom  = 0.5 * (pre.momentum + post.momentum);
0624         double   path_len = (post.position - pre.position).R();
0625         Geant4Tracker::Hit* hit = new Geant4Tracker::Hit(pre.truth.trackID, pre.truth.pdgID,
0626                                                          depo, time, path_len, pos, mom);
0627         hit->cellID   = cell;
0628         collection->add(hit);
0629         sensitive->printM2("+++ TrackID:%6d [%s] CREATE hit combination with %2d deposit(s):"
0630                            " %e MeV  Pos:%8.2f %8.2f %8.2f",
0631                            pre.truth.trackID,sensitive->c_name(),combined,pre.truth.deposit/CLHEP::MeV,
0632                            pos.X()/CLHEP::mm,pos.Y()/CLHEP::mm,pos.Z()/CLHEP::mm);
0633         clear();
0634       }
0635 
0636 
0637       /// Method for generating hit(s) using the information of G4Step object.
0638       G4bool process(const G4Step* step, G4TouchableHistory* ) {
0639         Geant4StepHandler h(step);
0640         // std::cout << " process called - pre pos: " << h.prePos() << " post pos " << h.postPos() 
0641         //    << " edep: " << h.deposit() << std::endl ;
0642         void *prePV = h.volume(h.pre), *postPV = h.volume(h.post);
0643 
0644         Geant4HitCollection* coll = sensitive->collection(0);
0645         /// If we are handling a new track, then store the content of the previous one.
0646         if ( mustSaveTrack(h.track) )  {
0647           extractHit(coll);
0648         }
0649         /// Initialize the deposits of the next hit.
0650         if ( current < 0 )  {
0651           start(step, h.pre);
0652         }
0653         /// ....update .....
0654         update(h);
0655 
0656         if ( prePV != postPV ) {
0657           void* postSD = h.sd(h.post);
0658           extractHit(coll);
0659           if ( 0 != postSD )   {
0660             void* preSD = h.sd(h.pre);
0661             if ( preSD == postSD ) {
0662               start(step,h.post);
0663             }
0664           }
0665         }
0666         else if ( h.track->GetTrackStatus() == fStopAndKill ) {
0667           extractHit(coll);
0668         }
0669         return true;
0670       }
0671 
0672       /// Method for generating hit(s) using the information of fast simulation spot object.
0673       G4bool process(const Geant4FastSimSpot* spot, G4TouchableHistory* ) {
0674         Geant4FastSimHandler h(spot);
0675         G4VPhysicalVolume*   prePV = firstSpotVolume, *postPV = h.volume();
0676         Geant4HitCollection* coll  = sensitive->collection(0);
0677         /// If we are handling a new track, then store the content of the previous one.
0678         if ( mustSaveTrack(h.track) )  {
0679           extractHit(coll);
0680         }
0681         /// Initialize the deposits of the next hit.
0682         if ( current < 0 )  {
0683           start(spot);
0684         }
0685         /// ....update .....
0686         update(h);
0687 
0688         if ( firstSpotVolume && prePV != postPV )   {
0689           void* postSD = h.sd();
0690           extractHit(coll);
0691           if ( 0 != postSD )   {
0692             void* preSD = prePV ? prePV->GetLogicalVolume()->GetSensitiveDetector() : nullptr;
0693             if ( preSD == postSD ) {
0694               start(spot);
0695             }
0696           }
0697         }
0698         else if ( h.track->GetTrackStatus() == fStopAndKill ) {
0699           extractHit(coll);
0700         }
0701         return true;
0702       }
0703 
0704       /// Post-event action callback
0705       void endEvent(const G4Event* /* event */)   {
0706         // We need to add the possibly last added hit to the collection here.
0707         // otherwise the last hit would be assigned to the next event and the
0708         // MC truth would be screwed.
0709         //
0710         // Alternatively the 'update' method would become rather CPU consuming,
0711         // beacuse the extract action would have to be recalculated over and over.
0712         if ( current > 0 )   {
0713           Geant4HitCollection* coll = sensitive->collection(0);
0714           extractHit(coll);
0715         }
0716       }
0717     };
0718 
0719     /// Initialization overload for specialization
0720     template <> void Geant4SensitiveAction<TrackerCombine>::initialize() {
0721       eventAction().callAtEnd(&m_userData,&TrackerCombine::endEvent);
0722       m_userData.e_cut = m_sensitive.energyCutoff();
0723       m_userData.sensitive = this;
0724     }
0725 
0726     /// Define collections created by this sensitivie action object
0727     template <> void Geant4SensitiveAction<TrackerCombine>::defineCollections() {
0728       m_collectionID = declareReadoutFilteredCollection<Geant4Tracker::Hit>();
0729     }
0730 
0731     /// Method for generating hit(s) using the information of G4Step object.
0732     template <> void Geant4SensitiveAction<TrackerCombine>::clear(G4HCofThisEvent*) {
0733       m_userData.clear();
0734     }
0735 
0736     /// Method for generating hit(s) using the information of G4Step object.
0737     template <> G4bool
0738     Geant4SensitiveAction<TrackerCombine>::process(const G4Step* step, G4TouchableHistory* history) {
0739       return m_userData.process(step, history);
0740     }
0741 
0742     /// GFlash/FastSim interface: Method for generating hit(s) using the information of Geant4FastSimSpot object.
0743     template <> bool
0744     Geant4SensitiveAction<TrackerCombine>::processFastSim(const Geant4FastSimSpot* spot,
0745                               G4TouchableHistory*      history)    {
0746       return m_userData.process(spot, history);
0747     }
0748 
0749     typedef Geant4SensitiveAction<TrackerCombine>  Geant4TrackerCombineAction;
0750     typedef Geant4TrackerAction                    Geant4SimpleTrackerAction;
0751     typedef Geant4CalorimeterAction                Geant4SimpleCalorimeterAction;
0752     typedef Geant4OpticalCalorimeterAction         Geant4SimpleOpticalCalorimeterAction;
0753   }
0754 }
0755 
0756 using namespace dd4hep::sim;
0757 
0758 #include <DDG4/Factories.h>
0759 // Special void entry point
0760 DECLARE_GEANT4SENSITIVE(Geant4VoidSensitiveAction)
0761 // Standard factories used for simulation
0762 DECLARE_GEANT4SENSITIVE(Geant4TrackerAction)
0763 DECLARE_GEANT4SENSITIVE(Geant4OpticalTrackerAction)
0764 DECLARE_GEANT4SENSITIVE(Geant4TrackerCombineAction)
0765 DECLARE_GEANT4SENSITIVE(Geant4CalorimeterAction)
0766 DECLARE_GEANT4SENSITIVE(Geant4OpticalCalorimeterAction)
0767 DECLARE_GEANT4SENSITIVE(Geant4ScintillatorCalorimeterAction)
0768 
0769 // Need these factories for backwards compatibility
0770 DECLARE_GEANT4SENSITIVE(Geant4SimpleTrackerAction)
0771 DECLARE_GEANT4SENSITIVE(Geant4SimpleCalorimeterAction)
0772 DECLARE_GEANT4SENSITIVE(Geant4SimpleOpticalCalorimeterAction)