Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
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
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
// $Id: Geant4Converter.cpp 603 2013-06-13 21:15:14Z markus.frank $
//====================================================================
// AIDA Detector description implementation for LCD
//--------------------------------------------------------------------
//
// Author : M.Frank
//
//====================================================================
// Framework include files
#include "DD4hep/Printout.h"
#include "DD4hep/Primitives.h"
#include "DD4hep/InstanceCount.h"
#include "DDG4/Geant4Kernel.h"
#include "DDG4/Geant4Mapping.h"
#include "DDG4/Geant4SensDetAction.h"
#include "DDG4/Geant4StepHandler.h"
#include "DDG4/Geant4VolumeManager.h"
#include <G4SDManager.hh>
#include <G4VSensitiveDetector.hh>
// C/C++ include files
#include <stdexcept>
using namespace std;
using namespace DD4hep;
using namespace DD4hep::Simulation;
namespace {
Geant4ActionSD* _getSensitiveDetector(const string& name) {
G4SDManager* mgr = G4SDManager::GetSDMpointer();
G4VSensitiveDetector* sd = mgr->FindSensitiveDetector(name);
if ( 0 == sd ) {
throw runtime_error(format("Geant4Sensitive","DDG4: You requested to configure actions "
"for the sensitive detector %s,\nDDG4: which is not known to Geant4. "
"Are you sure you already converted the geometry?",name.c_str()));
}
Geant4ActionSD* action_sd = dynamic_cast<Geant4ActionSD*>(sd);
if ( 0 == action_sd ) {
throw runtime_error(format("Geant4Sensitive","DDG4: You may only configure actions "
"for sensitive detectors of type Geant4ActionSD.\n"
"DDG4: The sensitive detector of %s is of type %s, which is incompatible.",
name.c_str(),typeinfoName(typeid(*sd)).c_str()));
}
return action_sd;
}
}
/// Standard action constructor
Geant4ActionSD::Geant4ActionSD(const std::string& name)
: Geant4Action(0,name)
{
InstanceCount::increment(this);
}
/// Default destructor
Geant4ActionSD::~Geant4ActionSD()
{
InstanceCount::decrement(this);
}
/// Standard constructor
Geant4Filter::Geant4Filter(Geant4Context* context, const std::string& name)
: Geant4Action(context,name)
{
InstanceCount::increment(this);
}
/// Standard destructor
Geant4Filter::~Geant4Filter() {
InstanceCount::decrement(this);
}
/// Filter action. Return true if hits should be processed
bool Geant4Filter::operator()(const G4Step*) const {
return true;
}
/// Constructor. The detector element is identified by the name
Geant4Sensitive::Geant4Sensitive(Geant4Context* ctxt, const string& name, DetElement det, LCDD& lcdd)
: Geant4Action(ctxt,name), m_sensitiveDetector(0), m_sequence(0), m_lcdd(lcdd),
m_detector(det), m_sensitive(), m_readout()
{
InstanceCount::increment(this);
if ( !det.isValid() ) {
throw runtime_error(format("Geant4Sensitive","DDG4: Detector elemnt for %s is invalid.",name.c_str()));
}
m_sequence = ctxt->kernel().sensitiveAction(m_detector.name());
m_sensitive = lcdd.sensitiveDetector(det.name());
m_readout = m_sensitive.readout();
}
/// Standard destructor
Geant4Sensitive::~Geant4Sensitive() {
m_filters(&Geant4Filter::release);
m_filters.clear();
InstanceCount::decrement(this);
}
/// Add an actor responding to all callbacks. Sequence takes ownership.
void Geant4Sensitive::adopt(Geant4Filter* filter) {
if ( filter ) {
filter->addRef();
m_filters.add(filter);
return;
}
throw runtime_error("Geant4SensDetActionSequence: Attempt to add invalid sensitive filter!");
}
/// Callback before hit processing starts. Invoke all filters.
bool Geant4Sensitive::accept(const G4Step* step) const {
bool result = m_filters.filter(&Geant4Filter::operator(),step);
return result;
}
/// Access to the sensitive detector object
void Geant4Sensitive::setDetector(Geant4ActionSD* sens_det) {
m_sensitiveDetector = sens_det;
}
/// Access to the sensitive detector object
Geant4ActionSD& Geant4Sensitive::detector() const {
if ( m_sensitiveDetector ) return *m_sensitiveDetector;
//m_sensitiveDetector = _getSensitiveDetector(m_detector.name());
//if ( m_sensitiveDetector ) return *m_sensitiveDetector;
throw runtime_error(format("Geant4Sensitive","DDG4: The sensitive detector for action %s "
"was not properly configured.",name().c_str()));
}
/// Access to the hosting sequence
Geant4SensDetActionSequence& Geant4Sensitive::sequence() const {
return *m_sequence;
}
/// Access HitCollection container names
const string& Geant4Sensitive::hitCollectionName(size_t which) const {
return sequence().hitCollectionName(which);
}
/// Retrieve the hits collection associated with this detector by its serial number
Geant4HitCollection* Geant4Sensitive::collection(size_t which) {
return sequence().collection(which);
}
/// Retrieve the hits collection associated with this detector by its collection identifier
Geant4HitCollection* Geant4Sensitive::collectionByID(size_t id) {
return sequence().collectionByID(id);
}
/// Method invoked at the begining of each event.
void Geant4Sensitive::begin(G4HCofThisEvent* /* HCE */) {
}
/// Method invoked at the end of each event.
void Geant4Sensitive::end(G4HCofThisEvent* /* HCE */) {
}
/// Method for generating hit(s) using the information of G4Step object.
bool Geant4Sensitive::process(G4Step* /* step */, G4TouchableHistory* /* history */) {
return false;
}
/// Method is invoked if the event abortion is occured.
void Geant4Sensitive::clear(G4HCofThisEvent* /* HCE */) {
}
/// Returns the volumeID of the sensitive volume corresponding to the step
long long int Geant4Sensitive::volumeID(G4Step* s) {
Geant4StepHandler step(s);
Geant4VolumeManager volMgr = Geant4Mapping::instance().volumeManager();
VolumeID id = volMgr.volumeID(step.preTouchable());
return id;
}
/// Standard constructor
Geant4SensDetActionSequence::Geant4SensDetActionSequence(Geant4Context* context,const string& nam)
: Geant4Action(context,nam), m_hce(0)
{
InstanceCount::increment(this);
context->sensitiveActions().insert(name(),this);
}
/// Default destructor
Geant4SensDetActionSequence::~Geant4SensDetActionSequence() {
m_filters(&Geant4Filter::release);
m_actors(&Geant4Sensitive::release);
m_filters.clear();
m_actors.clear();
InstanceCount::decrement(this);
}
/// Add an actor responding to all callbacks
void Geant4SensDetActionSequence::adopt(Geant4Sensitive* sensitive) {
if ( sensitive ) {
sensitive->addRef();
m_actors.add(sensitive);
return;
}
throw runtime_error("Geant4SensDetActionSequence: Attempt to add invalid sensitive actor!");
}
/// Add an actor responding to all callbacks. Sequence takes ownership.
void Geant4SensDetActionSequence::adopt(Geant4Filter* filter) {
if ( filter ) {
filter->addRef();
m_filters.add(filter);
return;
}
throw runtime_error("Geant4SensDetActionSequence: Attempt to add invalid sensitive filter!");
}
/// Initialize the usage of a hit collection. Returns the collection identifier
size_t Geant4SensDetActionSequence::defineCollection(const std::string& name,create_t func) {
m_collections.push_back(make_pair(name,func));
return m_collections.size()-1;
}
/// Called at construction time of the sensitive detector to declare all hit collections
size_t Geant4SensDetActionSequence::Geant4SensDetActionSequence::defineCollections(Geant4ActionSD* sens_det) {
size_t count = 0;
m_detector = sens_det;
m_actors(&Geant4Sensitive::setDetector,sens_det);
for(HitCollections::const_iterator i=m_collections.begin(); i!=m_collections.end(); ++i) {
sens_det->defineCollection((*i).first);
++count;
}
return count;
}
/// Access HitCollection container names
const std::string& Geant4SensDetActionSequence::hitCollectionName(size_t which) const {
if ( which < m_collections.size() ) {
return m_collections[which].first;
}
static string blank="";
except("The collection name index for subdetector %s is out of range!",c_name());
return blank;
}
/// Retrieve the hits collection associated with this detector by its serial number
Geant4HitCollection* Geant4SensDetActionSequence::collection(size_t which) const {
if ( which < m_collections.size() ) {
int hc_id = m_detector->GetCollectionID(which);
Geant4HitCollection* c = (Geant4HitCollection*)m_hce->GetHC(hc_id);
if ( c ) return c;
except("The collection index for subdetector %s is wrong!",c_name());
}
except("The collection name index for subdetector %s is out of range!",c_name());
return 0;
}
/// Retrieve the hits collection associated with this detector by its collection identifier
Geant4HitCollection* Geant4SensDetActionSequence::collectionByID(size_t id) const {
Geant4HitCollection* c = (Geant4HitCollection*)m_hce->GetHC(id);
if ( c ) return c;
except("The collection index for subdetector %s is wrong!",c_name());
return 0;
}
/// Callback before hit processing starts. Invoke all filters.
bool Geant4SensDetActionSequence::accept(const G4Step* step) const {
bool result = m_filters.filter(&Geant4Filter::operator(),step);
return result;
}
/// Function to process hits
bool Geant4SensDetActionSequence::process(G4Step* step, G4TouchableHistory* hist) {
bool result = false;
for(vector<Geant4Sensitive*>::iterator i=m_actors->begin(); i != m_actors->end(); ++i) {
Geant4Sensitive* s = *i;
if ( s->accept(step) ) result |= s->process(step,hist);
}
m_process(step,hist);
return result;
}
/** G4VSensitiveDetector interface: Method invoked at the begining of each event.
* The hits collection(s) created by this sensitive detector must
* be set to the G4HCofThisEvent object at one of these two methods.
*/
void Geant4SensDetActionSequence::begin(G4HCofThisEvent* hce) {
m_hce = hce;
for(size_t count = 0; count < m_collections.size(); ++count) {
const std::pair<string,create_t>& cr = m_collections[count];
G4VHitsCollection* c = (*cr.second)(name(), cr.first);
int id = m_detector->GetCollectionID(count);
m_hce->AddHitsCollection(id,c);
}
m_actors(&Geant4Sensitive::begin,m_hce);
m_begin(m_hce);
}
/// G4VSensitiveDetector interface: Method invoked at the end of each event.
void Geant4SensDetActionSequence::end(G4HCofThisEvent* hce) {
m_end(hce);
m_actors(&Geant4Sensitive::end,hce);
m_hce = 0;
}
/// G4VSensitiveDetector interface: Method invoked if the event was aborted.
/** Hits collections created but not beibg set to G4HCofThisEvent
* at the event should be deleted.
* Collection(s) which have already set to G4HCofThisEvent
* will be deleted automatically.
*/
void Geant4SensDetActionSequence::clear() {
m_clear(m_hce);
m_actors(&Geant4Sensitive::clear,m_hce);
}
/// Default constructor
Geant4SensDetSequences::Geant4SensDetSequences() {
}
/// Default destructor
Geant4SensDetSequences::~Geant4SensDetSequences() {
for_each(m_sequences.begin(),m_sequences.end(),releaseObjects(m_sequences));
m_sequences.clear();
}
/// Access sequence member by name
Geant4SensDetActionSequence* Geant4SensDetSequences::operator[](const string& name) const {
string nam = "SD_Seq_"+name;
Members::const_iterator i=m_sequences.find(nam);
if ( i != m_sequences.end() )
return (*i).second;
throw runtime_error("Attempt to access undefined SensDetActionSequence!");
}
/// Access sequence member by name
Geant4SensDetActionSequence* Geant4SensDetSequences::find(const std::string& name) const {
string nam = "SD_Seq_"+name;
Members::const_iterator i=m_sequences.find(nam);
if ( i != m_sequences.end() )
return (*i).second;
return 0;
}
/// Insert sequence member
void Geant4SensDetSequences::insert(const string& name,Geant4SensDetActionSequence* seq) {
if ( seq ) {
string nam = "SD_Seq_"+name;
seq->addRef();
m_sequences[nam] = seq;
return;
}
throw runtime_error(format("Geant4SensDetSequences","Attempt to add invalid sensitive "
"sequence with name:%s",name.c_str()));
}