Newer
Older
Markus Frank
committed
// $Id: Handle.h 570 2013-05-17 07:47:11Z markus.frank $
//==========================================================================
// AIDA Detector description implementation for LCD
Markus Frank
committed
//--------------------------------------------------------------------------
// Copyright (C) Organisation européenne pour la Recherche nucléaire (CERN)
// All rights reserved.
Markus Frank
committed
// For the licensing terms see $DD4hepINSTALL/LICENSE.
// For the list of contributors see $DD4hepINSTALL/doc/CREDITS.
Markus Frank
committed
// Author : M.Frank
//
//==========================================================================
#ifndef DD4HEP_DDG4_GEANT4ESCAPECOUNTER_H
#define DD4HEP_DDG4_GEANT4ESCAPECOUNTER_H
Markus Frank
committed
// Framework include files
#include "DD4hep/Detector.h"
#include "DDG4/Geant4SensDetAction.h"
#include "DDG4/Geant4SteppingAction.h"
/// Namespace for the AIDA detector description toolkit
/// Namespace for the Geant4 based simulation part of the AIDA detector description toolkit
/// Class to measure the energy of escaping tracks
Markus Frank
committed
/** Class to measure the energy of escaping tracks of a detector using Geant 4
* \author M.Frank
* \version 1.0
* \ingroup DD4HEP_SIMULATION
*/
class Geant4EscapeCounter : /* virtual public Geant4SteppingAction, virtual */ public Geant4Sensitive {
/// Collection identifiers
size_t m_collectionID;
std::vector<std::string> m_detectorNames;
public:
/// Standard constructor
Markus Frank
committed
Geant4EscapeCounter(Geant4Context* ctxt, const std::string& name, DetElement det, LCDD& lcdd);
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
/// Default destructor
virtual ~Geant4EscapeCounter();
/// G4VSensitiveDetector interface: Method for generating hit(s) using the information of G4Step object.
virtual bool process(G4Step* step, G4TouchableHistory* history);
};
} // End namespace Simulation
} // End namespace DD4hep
#endif /* DD4HEP_DDG4_GEANT4ESCAPECOUNTER_H */
// $Id: Geant4Converter.cpp 603 2013-06-13 21:15:14Z markus.frank $
//====================================================================
// AIDA Detector description implementation for LCD
//--------------------------------------------------------------------
//
// Author : M.Frank
//
//====================================================================
#include "DD4hep/Printout.h"
#include "DD4hep/InstanceCount.h"
#include "DDG4/Geant4TouchableHandler.h"
#include "DDG4/Geant4TrackHandler.h"
#include "DDG4/Geant4StepHandler.h"
#include "DDG4/Geant4Mapping.h"
#include "DDG4/Geant4Data.h"
#include "CLHEP/Units/SystemOfUnits.h"
#include "G4VProcess.hh"
using namespace std;
using namespace DD4hep;
using namespace DD4hep::Simulation;
/// Standard constructor
Markus Frank
committed
Geant4EscapeCounter::Geant4EscapeCounter(Geant4Context* ctxt, const string& nam, DetElement det, LCDD& lcdd_ref)
Markus Frank
committed
: Geant4Sensitive(ctxt, nam, det, lcdd_ref)
{
string coll_name = name()+"Hits";
m_needsControl = true;
declareProperty("Shells",m_detectorNames);
m_collectionID = defineCollection<SimpleTracker::Hit>(coll_name);
InstanceCount::increment(this);
}
/// Default destructor
Geant4EscapeCounter::~Geant4EscapeCounter() {
InstanceCount::decrement(this);
}
/// G4VSensitiveDetector interface: Method for generating hit(s) using the information of G4Step object.
bool Geant4EscapeCounter::process(G4Step* step, G4TouchableHistory* /* history */) {
Geant4TrackHandler th(h.track);
Geant4TouchableHandler handler(step);
Markus Frank
committed
string hdlr_path = handler.path();
Position prePos = h.prePos();
HitCollection* coll = collection(m_collectionID);
SimpleTracker::Hit* hit = new SimpleTracker::Hit(th.id(),th.pdgID(),h.deposit(),th.time());
hit->cellID = volumeID(step);
hit->energyDeposit = th.energy();
hit->position = prePos;
hit->momentum = h.trkMom();
hit->length = 0;
coll->add(hit);
mark(h.track);
Markus Frank
committed
print("+++ Track:%4d %8.2f MeV [%s] %s Geant4 path:%s",
Markus Frank
committed
h.trkID(),h.trkEnergy()/CLHEP::MeV,th.name().c_str(),
Markus Frank
committed
th.creatorName().c_str(),hdlr_path.c_str());
// Kill track, so that it does no longer participate in the propagation
h.track->SetTrackStatus(fStopAndKill);
return true;
}
#include "DDG4/Factories.h"
DECLARE_GEANT4SENSITIVE(Geant4EscapeCounter)