source: trunk/Documentation/geant4/UserDocumentation/UsersGuides/ForApplicationDeveloper/html/ch04s04.html@ 901

Last change on this file since 901 was 901, checked in by garnier, 17 years ago

Add Geant4 Documentation at 8.12.2008

File size: 44.3 KB
Line 
1<html><head><meta http-equiv="Content-Type" content="text/html; charset=ISO-8859-1"><title>4.4.  Hits</title><link rel="stylesheet" href="../xml/XSLCustomizationLayer/G4HTMLStylesheet.css" type="text/css"><meta name="generator" content="DocBook XSL Stylesheets V1.71.1"><link rel="start" href="index.html" title="Geant4 User's Guide for Application Developers"><link rel="up" href="ch04.html" title="Chapter 4.  Detector Definition and Response"><link rel="prev" href="ch04s03.html" title="4.3.  Electromagnetic Field"><link rel="next" href="ch04s05.html" title="4.5.  Digitization"><script language="JavaScript">
2function remote_win(fName)
3{
4 var url = "AllResources/Detector/geometry.src/" + fName;
5 RemoteWin=window.open(url,"","resizable=no,toolbar=0,location=0,directories=0,status=0,menubar=0,scrollbars=0,copyhistory=0,width=520,height=520")
6 RemoteWin.creator=self
7}
8</script></head><body bgcolor="white" text="black" link="#0000FF" vlink="#840084" alink="#0000FF"><div class="navheader"><table width="100%" summary="Navigation header"><tr><th colspan="3" align="center">4.4. 
9Hits
10</th></tr><tr><td width="20%" align="left"><a accesskey="p" href="ch04s03.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><th width="60%" align="center">Chapter 4. 
11Detector Definition and Response
12</th><td width="20%" align="right"> <a accesskey="n" href="ch04s05.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr></table><hr></div><div class="sect1" lang="en"><div class="titlepage"><div><div><h2 class="title" style="clear: both"><a name="sect.Hits"></a>4.4. 
13Hits
14</h2></div></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.Hit"></a>4.4.1. 
15Hit
16</h3></div></div></div><p>
17A hit is a snapshot of the physical interaction of a track in the
18sensitive region of a detector. In it you can store information
19associated with a <span class="emphasis"><em>G4Step</em></span> object. This information can be
20
21</p><div class="itemizedlist"><ul type="disc" compact><li><p>
22 the position and time of the step,
23 </p></li><li><p>
24 the momentum and energy of the track,
25 </p></li><li><p>
26 the energy deposition of the step,
27 </p></li><li><p>
28 geometrical information,
29 </p></li></ul></div><p>
30
31or any combination of the above.
32</p><h5><a name="id427822"></a>
33G4VHit
34</h5><p>
35<span class="emphasis"><em>G4VHit</em></span> is an abstract base class which represents a hit.
36You must inherit this base class and derive your own concrete hit
37class(es). The member data of your concrete hit class can be, and
38should be, your choice.
39</p><p>
40<span class="emphasis"><em>G4VHit</em></span> has two virtual methods, <code class="literal">Draw()</code> and
41<code class="literal">Print()</code>. To draw or print out your concrete hits, these
42methods should be implemented. How to define the drawing method is
43described in <a href="ch08s09.html" title="8.9. 
44Polylines, Markers and Text
45">Section 8.9</a>.
46</p><h5><a name="id427865"></a>
47G4THitsCollection
48</h5><p>
49<span class="emphasis"><em>G4VHit</em></span> is an abstract class from which you derive your
50own concrete classes. During the processing of a given event,
51represented by a <span class="emphasis"><em>G4Event</em></span> object, many objects of the hit
52class will be produced, collected and associated with the event.
53Therefore, for each concrete hit class you must also prepare a
54concrete class derived from <span class="emphasis"><em>G4VHitsCollection</em></span>, an abstract
55class which represents a vector collection of user defined
56hits.
57</p><p>
58<span class="emphasis"><em>G4THitsCollection</em></span> is a template class derived from
59<span class="emphasis"><em>G4VHitsCollection</em></span>, and the concrete hit collection class of
60a particular <span class="emphasis"><em>G4VHit</em></span> concrete class can be instantiated from
61this template class. Each object of a hit collection must have a
62unique name for each event.
63</p><p>
64<span class="emphasis"><em>G4Event</em></span> has a <span class="emphasis"><em>G4HCofThisEvent</em></span> class
65object, that is a container class of collections of hits. Hit collections are
66stored by their pointers, whose type is that of the base class.
67</p><h5><a name="id427918"></a>
68An example of a concrete hit class
69</h5><p>
70<a href="ch04s04.html#programlist_Hits_1" title="Example 4.14. 
71An example of a concrete hit class.
72">Example 4.14</a> shows an example of a concrete hit class.
73
74</p><div class="example"><a name="programlist_Hits_1"></a><p class="title"><b>Example 4.14. 
75An example of a concrete hit class.
76</b></p><div class="example-contents"><pre class="programlisting">
77#ifndef ExN04TrackerHit_h
78#define ExN04TrackerHit_h 1
79
80#include "G4VHit.hh"
81#include "G4THitsCollection.hh"
82#include "G4Allocator.hh"
83#include "G4ThreeVector.hh"
84
85class ExN04TrackerHit : public G4VHit
86{
87 public:
88
89 ExN04TrackerHit();
90 ~ExN04TrackerHit();
91 ExN04TrackerHit(const ExN04TrackerHit &amp;right);
92 const ExN04TrackerHit&amp; operator=(const ExN04TrackerHit &amp;right);
93 int operator==(const ExN04TrackerHit &amp;right) const;
94
95 inline void * operator new(size_t);
96 inline void operator delete(void *aHit);
97
98 void Draw() const;
99 void Print() const;
100
101 private:
102 G4double edep;
103 G4ThreeVector pos;
104
105 public:
106 inline void SetEdep(G4double de)
107 { edep = de; }
108 inline G4double GetEdep() const
109 { return edep; }
110 inline void SetPos(G4ThreeVector xyz)
111 { pos = xyz; }
112 inline G4ThreeVector GetPos() const
113 { return pos; }
114
115};
116
117typedef G4THitsCollection&lt;ExN04TrackerHit&gt; ExN04TrackerHitsCollection;
118
119extern G4Allocator&lt;ExN04TrackerHit&gt; ExN04TrackerHitAllocator;
120
121inline void* ExN04TrackerHit::operator new(size_t)
122{
123 void *aHit;
124 aHit = (void *) ExN04TrackerHitAllocator.MallocSingle();
125 return aHit;
126}
127
128inline void ExN04TrackerHit::operator delete(void *aHit)
129{
130 ExN04TrackerHitAllocator.FreeSingle((ExN04TrackerHit*) aHit);
131}
132
133#endif
134</pre></div></div><p><br class="example-break">
135</p><p>
136<span class="emphasis"><em>G4Allocator</em></span> is a class for fast allocation of objects to
137the heap through the paging mechanism. For details of
138<span class="emphasis"><em>G4Allocator</em></span>, refer to <a href="ch03s02.html#sect.GeneManage" title="3.2.4. 
139General management classes
140">Section 3.2.4</a>.
141Use of <span class="emphasis"><em>G4Allocator</em></span>
142is not mandatory, but it is recommended, especially for users who
143are not familiar with the C++ memory allocation mechanism or
144alternative tools of memory allocation. On the other hand, note
145that <span class="emphasis"><em>G4Allocator</em></span> is to be used
146<span class="bold"><strong>only</strong></span> for the concrete
147class that is <span class="bold"><strong>not</strong></span> used as a base
148class of any other classes.
149For example, do <span class="bold"><strong>not</strong></span> use the
150<span class="emphasis"><em>G4Trajectory</em></span> class as a
151base class for a customized trajectory class, since
152<span class="emphasis"><em>G4Trajectory</em></span> uses <span class="emphasis"><em>G4Allocator</em></span>.
153</p><h5><a name="id428028"></a>
154G4THitsMap
155</h5><p>
156<span class="emphasis"><em>G4THitsMap</em></span> is an alternative to
157<span class="emphasis"><em>G4THitsCollection</em></span>.
158<span class="emphasis"><em>G4THitsMap</em></span> does not demand <span class="emphasis"><em>G4VHit</em></span>,
159but instead any variable which can be mapped with an integer key. Typically the key
160is a copy number of the volume, and the mapped value could for
161example be a double, such as the energy deposition in a volume.
162<span class="emphasis"><em>G4THitsMap</em></span> is convenient for applications which do not need
163to output event-by-event data but instead just accumulate them. All
164the <span class="emphasis"><em>G4VPrimitiveScorer</em></span> classes discussed in
165<a href="ch04s04.html#sect.Hits.G4Multi" title="4.4.5. 
166G4MultiFunctionalDetector and
167G4VPrimitiveScorer
168">Section 4.4.5</a> use <span class="emphasis"><em>G4THitsMap</em></span>.
169</p><p>
170<span class="emphasis"><em>G4THitsMap</em></span> is derived from the
171<span class="emphasis"><em>G4VHitsCollection</em></span>
172abstract base class and all objects of this class are also stored
173in <span class="emphasis"><em>G4HCofThisEvent</em></span> at the end of an event. How to access a
174<span class="emphasis"><em>G4THitsMap</em></span> object is discussed in the
175following section (<a href="ch04s04.html#sect.Hits.G4Multi" title="4.4.5. 
176G4MultiFunctionalDetector and
177G4VPrimitiveScorer
178">Section 4.4.5</a>).
179</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.SensDet"></a>4.4.2. 
180Sensitive detector
181</h3></div></div></div><h5><a name="id428112"></a>
182G4VSensitiveDetector
183</h5><p>
184<span class="emphasis"><em>G4VSensitiveDetector</em></span> is an abstract base class which
185represents a detector. The principal mandate of a sensitive
186detector is the construction of hit objects using information from
187steps along a particle track. The <code class="literal">ProcessHits()</code> method of
188<span class="emphasis"><em>G4VSensitiveDetector</em></span> performs this task using
189<span class="emphasis"><em>G4Step</em></span>
190objects as input. In the case of a "Readout" geometry (see
191<a href="ch04s04.html#sect.Hits.ReadGeom" title="4.4.3. 
192Readout geometry
193">Section 4.4.3</a>), objects of the
194<span class="emphasis"><em>G4TouchableHistory</em></span> class may be used as an optional input.
195</p><p>
196Your concrete detector class should be instantiated with the
197unique name of your detector. The name can be associated with one
198or more global names with "/" as a delimiter for categorizing your
199detectors. For example
200
201</p><div class="informalexample"><pre class="programlisting">
202 myEMcal = new MyEMcal("/myDet/myCal/myEMcal");
203</pre></div><p>
204
205where <code class="literal">myEMcal</code> is the name of your detector. The pointer to
206your sensitive detector must be set to one or more
207<span class="emphasis"><em>G4LogicalVolume</em></span> objects to set the sensitivity of these
208volumes. The pointer should also be registered to
209<span class="emphasis"><em>G4SDManager</em></span>, as described in
210<a href="ch04s04.html#sect.Hits.G4SDMan" title="4.4.4. 
211G4SDManager
212">Section 4.4.4</a>.
213</p><p>
214<span class="emphasis"><em>G4VSensitiveDetector</em></span> has three major virtual methods.
215
216</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
217 <code class="literal">ProcessHits()</code>
218 </span></dt><dd><p>
219 This method is invoked by <span class="emphasis"><em>G4SteppingManager</em></span> when a step
220 is composed in the <span class="emphasis"><em>G4LogicalVolume</em></span> which has the pointer
221 to this sensitive detector. The first argument of this method is a
222 <span class="emphasis"><em>G4Step</em></span> object of the current step. The second argument is a
223 <span class="emphasis"><em>G4TouchableHistory</em></span> object for the ``Readout geometry''
224 described in the next section. The second argument is <code class="literal">NULL</code>
225 if ``Readout geometry'' is not assigned to this sensitive detector.
226 In this method, one or more <span class="emphasis"><em>G4VHit</em></span> objects should be
227 constructed if the current step is meaningful for your
228 detector.
229 </p></dd><dt><span class="term">
230 <code class="literal">Initialize()</code>
231 </span></dt><dd><p>
232 This method is invoked at the beginning of each event. The
233 argument of this method is an object of the <span class="emphasis"><em>G4HCofThisEvent</em></span>
234 class. Hit collections, where hits produced in this particular
235 event are stored, can be associated with the <span class="emphasis"><em>G4HCofThisEvent</em></span>
236 object in this method. The hit collections associated with the
237 <span class="emphasis"><em>G4HCofThisEvent</em></span> object during this method can be used for
238 ``during the event processing'' digitization.
239 </p></dd><dt><span class="term">
240 <code class="literal">EndOfEvent()</code>
241 </span></dt><dd><p>
242 This method is invoked at the end of each event. The argument
243 of this method is the same object as the previous method. Hit
244 collections occasionally created in your sensitive detector can be
245 associated with the <span class="emphasis"><em>G4HCofThisEvent</em></span> object.
246 </p></dd></dl></div><p>
247</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.ReadGeom"></a>4.4.3. 
248Readout geometry
249</h3></div></div></div><p>
250This section describes how a ``Readout geometry'' can be defined. A
251Readout geometry is a virtual, parallel geometry for obtaining the
252channel number.
253</p><p>
254As an example, the accordion calorimeter of <span class="bold"><strong>ATLAS</strong></span>
255has a complicated tracking geometry, however the readout can be done by
256simple cylindrical sectors divided by theta, phi, and depth. Tracks
257will be traced in the tracking geometry, the ``real'' one, and the
258sensitive detector will have its own readout geometry Geant4 will
259message to find to which ``readout'' cell the current hit
260belongs.
261
262</p><div class="figure"><a name="fig.Hits_1"></a><div class="figure-contents"><div class="mediaobject" align="center"><img src="./AllResources/Detector/hit.src/RO.gif" align="middle" alt="Association of tracking and readout geometry."></div></div><p class="title"><b>Figure 4.8. 
263Association of tracking and readout geometry.
264</b></p></div><p><br class="figure-break">
265</p><p>
266<a href="ch04s04.html#fig.Hits_1" title="Figure 4.8. 
267Association of tracking and readout geometry.
268">Figure 4.8</a> shows how this association is done in Geant4.
269The first step is to associate a sensitive detector to a volume of the
270tracking geometry, in the usual way (see
271<a href="ch04s04.html#sect.Hits.SensDet" title="4.4.2. 
272Sensitive detector
273">Section 4.4.2</a>). The next step is to associate your
274<span class="emphasis"><em>G4VReadoutGeometry</em></span> object to the sensitive detector.
275</p><p>
276At tracking time, the base class <span class="emphasis"><em>G4VReadoutGeometry</em></span> will
277provide to your sensitive detector code the
278<span class="emphasis"><em>G4TouchableHistory</em></span> in the Readout geometry at the beginning
279of the step position (position of <span class="emphasis"><em>PreStepPoint</em></span> of
280<span class="emphasis"><em>G4Step</em></span>) and at this position only.
281</p><p>
282This <span class="emphasis"><em>G4TouchableHistory</em></span> is given to your sensitive
283detector code through the <span class="emphasis"><em>G4VSensitiveDetector</em></span> virtual
284method:
285
286</p><div class="informalexample"><pre class="programlisting">
287 G4bool processHits(G4Step* aStep, G4TouchableHistory* ROhist);
288</pre></div><p>
289
290by the <code class="literal">ROhist</code> argument.
291</p><p>
292Thus, you will be able to use information from both the
293<span class="emphasis"><em>G4Step</em></span> and the <span class="emphasis"><em>G4TouchableHistory</em></span>
294coming from your
295Readout geometry. Note that since the association is done through a
296sensitive detector object, it is perfectly possible to have several
297Readout geometries in parallel.
298</p><h5><a name="id428454"></a>
299Definition of a virtual geometry setup
300</h5><p>
301The base class for the implementation of a Readout geometry is
302<span class="emphasis"><em>G4VReadoutGeometry</em></span>. This class has a single pure virtual
303protected method:
304
305</p><div class="informalexample"><pre class="programlisting">
306 virtual G4VPhysicalVolume* build() = 0;
307</pre></div><p>
308
309which you must override in your concrete class. The
310<span class="emphasis"><em>G4VPhysicalVolume</em></span> pointer you will have to return is of the
311physical world of the Readout geometry.
312</p><p>
313The step by step procedure for constructing a Readout geometry is:
314
315</p><div class="itemizedlist"><ul type="disc" compact><li><p>
316 inherit from <span class="emphasis"><em>G4VReadoutGeometry</em></span> to define a
317 <span class="emphasis"><em>MyROGeom</em></span> class;
318 </p></li><li><p>
319 implement the Readout geometry in the <code class="literal">build()</code> method,
320 returning the physical world of this geometry.
321 </p><p>
322 The world is specified in the same way as for the detector
323 construction: a physical volume with no mother. The axis system of
324 this world is the same as the one of the world for tracking.
325 </p><p>
326 </p><p>
327 In this geometry you need to declare the sensitive parts in the
328 same way as in the tracking geometry: by setting a
329 non-<code class="literal">NULL</code> <span class="emphasis"><em>G4VSensitiveDetector</em></span>
330 pointer in, say, the
331 relevant <span class="emphasis"><em>G4LogicalVolume</em></span> objects. This sensitive class needs
332 to be there, but will not be used.
333 </p><p>
334 </p><p>
335 You will also need to assign well defined materials to the
336 volumes you place in this geometry, but these materials are
337 irrelevant since they will not be seen by the tracking. It is
338 foreseen to allow the setting of a <code class="literal">NULL</code> pointer in this
339 case of the parallel geometry.
340 </p><p>
341 </p></li><li><p>
342 in the <code class="literal">construct()</code> method of your concrete
343 <span class="emphasis"><em>G4VUserDetectorConstruction</em></span> class:
344 </p><p>
345 </p><div class="itemizedlist"><ul type="circle" compact><li><p>
346 instantiate your Readout geometry:
347 </p><div class="informalexample"><pre class="programlisting">
348 MyROGeom* ROgeom = new MyROGeom("ROName");
349 </pre></div><p>
350 </p></li><li><p>
351 build it:
352 </p><div class="informalexample"><pre class="programlisting">
353 ROgeom-&gt;buildROGeometry();
354 </pre></div><p>
355 That will invoke your <code class="literal">build()</code> method.
356 </p></li><li><p>
357 Instantiate the sensitive detector which will receive the
358 <code class="literal">ROGeom</code> pointer, <code class="literal">MySensitive</code>,
359 and add this sensitive to the <span class="emphasis"><em>G4SDManager</em></span>.
360 Associate this sensitive to
361 the volume(s) of the tracking geometry as usual.
362 </p></li><li><p>
363 Associate the sensitive to the Readout geometry:
364 </p><div class="informalexample"><pre class="programlisting">
365 MySensitive-&gt;SetROgeometry(ROgeom);
366 </pre></div><p>
367 </p></li></ul></div><p>
368 </p><p>
369 </p></li></ul></div><p>
370</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.G4SDMan"></a>4.4.4. 
371G4SDManager
372</h3></div></div></div><p>
373<span class="emphasis"><em>G4SDManager</em></span> is the singleton manager class for sensitive
374detectors.
375</p><h5><a name="id428671"></a>
376Activation/inactivation of sensitive detectors
377</h5><p>
378The user interface commands <code class="literal">activate</code> and
379<code class="literal">inactivate</code> are available to control your sensitive
380detectors. For example:
381
382</p><div class="informalexample"><pre class="programlisting">
383/hits/activate detector_name
384/hits/inactivate detector_name
385</pre></div><p>
386
387where <code class="literal">detector_name</code> can be the detector name or the
388category name.
389</p><p>
390For example, if your EM calorimeter is named
391
392</p><div class="informalexample"><pre class="programlisting">
393/myDet/myCal/myEMcal
394/hits/inactivate myCal
395</pre></div><p>
396
397will inactivate all detectors belonging to the <code class="literal">myCal</code>
398category.
399</p><h5><a name="id428731"></a>
400Access to the hit collections
401</h5><p>Hit collections are accessed for various cases.
402
403</p><div class="itemizedlist"><ul type="disc" compact><li><p>
404 Digitization
405 </p></li><li><p>
406 Event filtering in <span class="emphasis"><em>G4VUserStackingAction</em></span>
407 </p></li><li><p>
408 ``End of event'' simple analysis
409 </p></li><li><p>
410 Drawing / printing hits
411 </p></li></ul></div><p>
412</p><p>
413The following is an example of how to access the hit collection
414of a particular concrete type:
415
416</p><div class="informalexample"><pre class="programlisting">
417 G4SDManager* fSDM = G4SDManager::GetSDMpointer();
418 G4RunManager* fRM = G4RunManager::GetRunManager();
419 G4int collectionID = fSDM-&gt;GetCollectionID("collection_name");
420 const G4Event* currentEvent = fRM-&gt;GetCurrentEvent();
421 G4HCofThisEvent* HCofEvent = currentEvent-&gt;GetHCofThisEvent();
422 <span class="emphasis"><em>MyHitsCollection</em></span>* myCollection = (<span class="emphasis"><em>MyHitsCollection</em></span>*)(HC0fEvent-&gt;GetHC(collectionID));
423</pre></div><p>
424</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.G4Multi"></a>4.4.5. 
425<span class="emphasis"><em>G4MultiFunctionalDetector</em></span> and
426<span class="emphasis"><em>G4VPrimitiveScorer</em></span>
427</h3></div></div></div><p>
428<span class="emphasis"><em>G4MultiFunctionalDetector</em></span> is a concrete class derived from
429<span class="emphasis"><em>G4VSensitiveDetector</em></span>. Instead of implementing a
430user-specific detector class, <span class="emphasis"><em>G4MultiFunctionalDetector</em></span>
431allows the user to register <span class="emphasis"><em>G4VPrimitiveScorer</em></span> classes to
432build up the sensitivity. <span class="emphasis"><em>G4MultiFunctionalDetector</em></span> should
433be instantiated in the users detector construction with its unique
434name and should be assigned to one or more <span class="emphasis"><em>G4LogicalVolume</em></span>s.
435</p><p>
436<span class="emphasis"><em>G4VPrimitiveScorer</em></span> is an abstract base class representing
437a class to be registered to <span class="emphasis"><em>G4MultiFunctionalDetector</em></span> that
438creates a <span class="emphasis"><em>G4THitsMap</em></span> object of one physics quantity for an
439event. Geant4 provides many concrete primitive scorer classes
440listed in <a href="ch04s04.html#sect.Hits.G4VPrim" title="4.4.6. 
441Concrete classes of G4VPrimitiveScorer
442">Section 4.4.6</a>, and the user can
443also implement his/her own primitive scorers. Each primitive scorer
444object must be instantiated with a name that must be unique among
445primitive scorers registered in a <span class="emphasis"><em>G4MultiFunctionalDetector</em></span>.
446Please note that a primitive scorer object must <span class="bold"><strong>not</strong></span> be
447shared by more than one <span class="emphasis"><em>G4MultiFunctionalDetector</em></span>
448object.
449</p><p>As mentioned in <a href="ch04s04.html#sect.Hits.Hit" title="4.4.1. 
450Hit
451">Section 4.4.1</a>,
452each <span class="emphasis"><em>G4VPrimitiveScorer</em></span> generates one
453<span class="emphasis"><em>G4THitsMap</em></span> object
454per event. The name of the map object is the same as the name of
455the primitive scorer. Each of the concrete primitive scorers listed
456in <a href="ch04s04.html#sect.Hits.G4VPrim" title="4.4.6. 
457Concrete classes of G4VPrimitiveScorer
458">Section 4.4.6</a> generates a
459<span class="emphasis"><em>G4THitsMap&lt;G4double&gt;</em></span> that maps a
460<span class="emphasis"><em>G4double</em></span> value
461to its key integer number. By default, the key is taken as the copy
462number of the <span class="emphasis"><em>G4LogicalVolume</em></span> to which
463<span class="emphasis"><em>G4MultiFunctionalDetector</em></span> is assigned. In case the logical
464volume is uniquely placed in its mother volume and the mother is
465replicated, the copy number of its mother volume can be taken by
466setting the second argument of the <span class="emphasis"><em>G4VPrimitiveScorer</em></span>
467constructor, "<span class="emphasis"><em>depth</em></span>" to 1, i.e. one level up. Furthermore,
468in case the key must consider more than one copy number of a
469different geometry hierarchy, the user can derive his/her own
470primitive scorer from the provided concrete class and implement the
471<span class="emphasis"><em>GetIndex(G4Step*)</em></span> virtual method to return the unique
472key.
473</p><p>
474<a href="ch04s04.html#programlist_Hits_2" title="Example 4.15. 
475An example of defining primitive sensitivity classes taken from
476ExN07DetectorConstruction.
477">Example 4.15</a> shows an example of primitive sensitivity
478class definitions.
479
480</p><div class="example"><a name="programlist_Hits_2"></a><p class="title"><b>Example 4.15. 
481An example of defining primitive sensitivity classes taken from
482<span class="emphasis"><em>ExN07DetectorConstruction</em></span>.
483</b></p><div class="example-contents"><pre class="programlisting">
484void ExN07DetectorConstruction::SetupDetectors()
485{
486 G4String filterName, particleName;
487
488 G4SDParticleFilter* gammaFilter =
489 new G4SDParticleFilter(filterName="gammaFilter",particleName="gamma");
490 G4SDParticleFilter* electronFilter =
491 new G4SDParticleFilter(filterName="electronFilter",particleName="e-");
492 G4SDParticleFilter* positronFilter =
493 new G4SDParticleFilter(filterName="positronFilter",particleName="e+");
494 G4SDParticleFilter* epFilter = new G4SDParticleFilter(filterName="epFilter");
495 epFilter-&gt;add(particleName="e-");
496 epFilter-&gt;add(particleName="e+");
497
498
499 for(G4int i=0;i&lt;3;i++)
500 {
501 for(G4int j=0;j&lt;2;j++)
502 {
503 // Loop counter j = 0 : absorber
504 // = 1 : gap
505 G4String detName = calName[i];
506 if(j==0)
507 { detName += "_abs"; }
508 else
509 { detName += "_gap"; }
510 G4MultiFunctionalDetector* det = new G4MultiFunctionalDetector(detName);
511
512 // The second argument in each primitive means the "level" of geometrical hierarchy,
513 // the copy number of that level is used as the key of the G4THitsMap.
514 // For absorber (j = 0), the copy number of its own physical volume is used.
515 // For gap (j = 1), the copy number of its mother physical volume is used, since there
516 // is only one physical volume of gap is placed with respect to its mother.
517 G4VPrimitiveScorer* primitive;
518 primitive = new G4PSEnergyDeposit("eDep",j);
519 det-&gt;RegisterPrimitive(primitive);
520 primitive = new G4PSNofSecondary("nGamma",j);
521 primitive-&gt;SetFilter(gammaFilter);
522 det-&gt;RegisterPrimitive(primitive);
523 primitive = new G4PSNofSecondary("nElectron",j);
524 primitive-&gt;SetFilter(electronFilter);
525 det-&gt;RegisterPrimitive(primitive);
526 primitive = new G4PSNofSecondary("nPositron",j);
527 primitive-&gt;SetFilter(positronFilter);
528 det-&gt;RegisterPrimitive(primitive);
529 primitive = new G4PSMinKinEAtGeneration("minEkinGamma",j);
530 primitive-&gt;SetFilter(gammaFilter);
531 det-&gt;RegisterPrimitive(primitive);
532 primitive = new G4PSMinKinEAtGeneration("minEkinElectron",j);
533 primitive-&gt;SetFilter(electronFilter);
534 det-&gt;RegisterPrimitive(primitive);
535 primitive = new G4PSMinKinEAtGeneration("minEkinPositron",j);
536 primitive-&gt;SetFilter(positronFilter);
537 det-&gt;RegisterPrimitive(primitive);
538 primitive = new G4PSTrackLength("trackLength",j);
539 primitive-&gt;SetFilter(epFilter);
540 det-&gt;RegisterPrimitive(primitive);
541 primitive = new G4PSNofStep("nStep",j);
542 primitive-&gt;SetFilter(epFilter);
543 det-&gt;RegisterPrimitive(primitive);
544
545 G4SDManager::GetSDMpointer()-&gt;AddNewDetector(det);
546 if(j==0)
547 { layerLogical[i]-&gt;SetSensitiveDetector(det); }
548 else
549 { gapLogical[i]-&gt;SetSensitiveDetector(det); }
550 }
551 }
552}
553</pre></div></div><p><br class="example-break">
554</p><p>
555Each <span class="emphasis"><em>G4THitsMap</em></span> object can be accessed from
556<span class="emphasis"><em>G4HCofThisEvent</em></span> with a unique collection ID number. This ID
557number can be obtained from <span class="emphasis"><em>G4SDManager::GetCollectionID()</em></span>
558with a name of <span class="emphasis"><em>G4MultiFunctionalDetector</em></span> and
559<span class="emphasis"><em>G4VPrimitiveScorer</em></span> connected with a slush ("/").
560<span class="emphasis"><em>G4THitsMap</em></span> has a [] operator taking the key value as an
561argument and returning <span class="bold"><strong>the pointer</strong></span> of the value.
562Please note that the [] operator returns
563<span class="bold"><strong>the pointer</strong></span> of the value. If
564you get zero from the [] operator, it does <span class="bold"><strong>not</strong></span> mean the
565value is zero, but that the provided key does not exist. The value
566itself is accessible with an astarisk ("*"). It is advised to check
567the validity of the returned pointer before accessing the value.
568<span class="emphasis"><em>G4THitsMap</em></span> also has a += operator in order to accumulate
569event data into run data. <a href="ch04s04.html#programlist_Hits_3" title="Example 4.16. 
570An example of accessing to G4THitsMap objects.
571">Example 4.16</a> shows the use of
572<span class="emphasis"><em>G4THitsMap</em></span>.
573
574</p><div class="example"><a name="programlist_Hits_3"></a><p class="title"><b>Example 4.16. 
575An example of accessing to <span class="emphasis"><em>G4THitsMap</em></span> objects.
576</b></p><div class="example-contents"><pre class="programlisting">
577#include "ExN07Run.hh"
578#include "G4Event.hh"
579#include "G4HCofThisEvent.hh"
580#include "G4SDManager.hh"
581
582ExN07Run::ExN07Run()
583{
584 G4String detName[6] = {"Calor-A_abs","Calor-A_gap","Calor-B_abs","Calor-B_gap",
585 "Calor-C_abs","Calor-C_gap"};
586 G4String primNameSum[6] = {"eDep","nGamma","nElectron","nPositron","trackLength","nStep"};
587 G4String primNameMin[3] = {"minEkinGamma","minEkinElectron","minEkinPositron"};
588
589 G4SDManager* SDMan = G4SDManager::GetSDMpointer();
590 G4String fullName;
591 for(size_t i=0;i&lt;6;i++)
592 {
593 for(size_t j=0;j&lt;6;j++)
594 {
595 fullName = detName[i]+"/"+primNameSum[j];
596 colIDSum[i][j] = SDMan-&gt;GetCollectionID(fullName);
597 }
598 for(size_t k=0;k&lt;3;k++)
599 {
600 fullName = detName[i]+"/"+primNameMin[k];
601 colIDMin[i][k] = SDMan-&gt;GetCollectionID(fullName);
602 }
603 }
604}
605
606
607void ExN07Run::RecordEvent(const G4Event* evt)
608{
609 G4HCofThisEvent* HCE = evt-&gt;GetHCofThisEvent();
610 if(!HCE) return;
611 numberOfEvent++;
612 for(size_t i=0;i&lt;6;i++)
613 {
614 for(size_t j=0;j&lt;6;j++)
615 {
616 G4THitsMap&lt;G4double&gt;* evtMap = (G4THitsMap&lt;G4double&gt;*)(HCE-&gt;GetHC(colIDSum[i][j]));
617 mapSum[i][j] += *evtMap;
618 }
619 for(size_t k=0;k&lt;3;k++)
620 {
621 G4THitsMap&lt;G4double&gt;* evtMap = (G4THitsMap&lt;G4double&gt;*)(HCE-&gt;GetHC(colIDMin[i][k]));
622 std::map&lt;G4int,G4double*&gt;::iterator itr = evtMap-&gt;GetMap()-&gt;begin();
623 for(; itr != evtMap-&gt;GetMap()-&gt;end(); itr++)
624 {
625 G4int key = (itr-&gt;first);
626 G4double val = *(itr-&gt;second);
627 G4double* mapP = mapMin[i][k][key];
628 if( mapP &amp;&amp; (val&gt;*mapP) ) continue;
629 mapMin[i][k].set(key,val);
630 }
631 }
632 }
633}
634</pre></div></div><p><br class="example-break">
635</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.G4VPrim"></a>4.4.6. 
636Concrete classes of <span class="emphasis"><em>G4VPrimitiveScorer</em></span>
637</h3></div></div></div><p>
638With Geant4 version 8.0, several concrete primitive scorer classes
639are provided, all of which are derived from the
640<span class="emphasis"><em>G4VPrimitiveScorer</em></span> abstract base class and which are to be
641registered to <span class="emphasis"><em>G4MultiFunctionalDetector</em></span>. Each of them
642contains one <span class="emphasis"><em>G4THitsMap</em></span> object and scores a simple double
643value for each key.
644</p><h5><a name="id429128"></a>
645Track length scorers
646</h5><p>
647</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
648 G4PSTrackLength
649 </span></dt><dd><p>
650 The track length is defined as the sum of step lengths of the
651 particles inside the cell. Bt default, the track weight is not
652 taken into account, but could be used as a multiplier of each step
653 length if the <span class="emphasis"><em>Weighted()</em></span> method of this class object is
654 invoked.
655 </p></dd><dt><span class="term">
656 G4PSPassageTrackLength
657 </span></dt><dd><p>
658 The passage track length is the same as the track length in
659 <span class="emphasis"><em>G4PSTrackLength</em></span>, except that only tracks which pass
660 through the volume are taken into account. It means newly-generated or
661 stopped tracks inside the cell are excluded from the calculation.
662 By default, the track weight is not taken into account, but could
663 be used as a multiplier of each step length if the
664 <span class="emphasis"><em>Weighted()</em></span> method of this class object is invoked.
665 </p></dd></dl></div><p>
666</p><h5><a name="id429183"></a>
667<span class="bold"><strong>Deposited energy scorers</strong></span>
668</h5><p>
669</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
670 G4PSEnergyDeposit
671 </span></dt><dd><p>
672 This scorer stores a sum of particles' energy deposits at each
673 step in the cell. The particle weight is multiplied at each
674 step.
675 </p></dd><dt><span class="term">
676 G4PSDoseDeposit
677 </span></dt><dd><p>
678 In some cases, dose is a more convenient way to evaluate the
679 effect of energy deposit in a cell than simple deposited energy.
680 The dose deposit is defined by the sum of energy deposits at each
681 step in a cell divided by the mass of the cell. The mass is
682 calculated from the density and volume of the cell taken from the
683 methods of <span class="emphasis"><em>G4VSolid</em></span> and
684 <span class="emphasis"><em>G4LogicalVolume</em></span>. The particle
685 weight is multiplied at each step.
686 </p></dd></dl></div><p>
687</p><h5><a name="id429239"></a>
688<span class="bold"><strong>Current and flux scorers</strong></span>
689</h5><p>
690There are two different definitions of a particle's flow for a
691given geometry. One is a current and the other is a flux. In our
692scorers, the current is simply defined as the number of particles
693(with the particle's weight) at a certain surface or volume, while
694the flux takes the particle's injection angle to the geometry into
695account. The current and flux are usually defined at a surface, but
696volume current and volume flux are also provided.
697</p><p>
698</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
699 G4PSFlatSurfaceCurrent
700 </span></dt><dd><p>
701 Flat surface current is a surface based scorer. The present
702 implementation is limited to scoring only at the -Z surface of a
703 <span class="emphasis"><em>G4Box</em></span> solid. The quantity is defined by the number
704 of tracks that reach the surface. The user must choose a direction of the
705 particle to be scored. The choices are fCurrent_In, fCurrent_Out,
706 or fCurrent_InOut, one of which must be entered as the second
707 argument of the constructor. Here, fCurrent_In scores incoming
708 particles to the cell, while fCurrent_Out scores only outgoing
709 particles from the cell. fCurrent_InOut scores both directions. The
710 current is multiplied by particle weight and is normalized for a
711 unit area.
712 </p></dd><dt><span class="term">
713 G4PSSphereSurfaceCurrent
714 </span></dt><dd><p>
715 Sphere surface current is a surface based scorer, and similar
716 to the G4PSFlatSurfaceCurrent. The only difference is that the
717 surface is defined at the <span class="bold"><strong>inner surface</strong></span>
718 of a G4Sphere solid.
719 </p></dd><dt><span class="term">
720 G4PSPassageCurrent
721 </span></dt><dd><p>
722 Passage current is a volume-based scorer. The current is
723 defined by the number of tracks that pass through the volume. A
724 particle weight is applied at the exit point. A passage current is
725 defined for a volume.
726 </p></dd><dt><span class="term">
727 G4PSFlatSurfaceFlux
728 </span></dt><dd><p>
729 Flat surface flux is a surface based flux scorer. The surface
730 flux is defined by the number of tracks that reach the surface. The
731 expression of surface flux is given by the sum of W/cos(t)/A, where
732 W, t and A represent particle weight, injection angle of particle
733 with respect to the surface normal, and area of the surface. The
734 user must enter one of the particle directions, fFlux_In,
735 fFlux_Out, or fFlux_InOut in the constructor. Here, fFlux_In scores
736 incoming particles to the cell, while fFlux_Out scores outgoing
737 particles from the cell. fFlux_InOut scores both directions.
738 </p></dd><dt><span class="term">
739 G4PSCellFlux
740 </span></dt><dd><p>
741 Cell flux is a volume based flux scorer. The cell flux is
742 defined by a track length (L) of the particle inside a volume
743 divided by the volume (V) of this cell. The track length is
744 calculated by a sum of the step lengths in the cell. The expression
745 for cell flux is given by the sum of (W*L)/V, where W is a particle
746 weight, and is multiplied by the track length at each step.
747 </p></dd><dt><span class="term">
748 G4PSPassageCellFlux
749 </span></dt><dd><p>
750 Passage cell flux is a volume based scorer similar to
751 <span class="emphasis"><em>G4PSCellFlux</em></span>. The only difference is that tracks
752 which pass through a cell are taken into account. It means generated or
753 stopped tracks inside the volume are excluded from the
754 calculation.
755 </p></dd></dl></div><p>
756</p><h5><a name="id429373"></a>
757<span class="bold"><strong>Other scorers</strong></span>
758</h5><p>
759</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
760 G4PSMinKinEAtGeneration
761 </span></dt><dd><p>
762 This scorer records the minimum kinetic energy of secondary
763 particles at their production point in the volume in an event. This
764 primitive scorer does not integrate the quantity, but records the
765 minimum quantity.
766 </p></dd><dt><span class="term">
767 G4PSNofSecondary
768 </span></dt><dd><p>
769 This class scores the number of secondary particles generated
770 in the volume. The weight of the secondary track is taken into
771 account.
772 </p></dd><dt><span class="term">
773 G4PSNofStep
774 </span></dt><dd><p>
775 This class scores the number of steps in the cell. A particle
776 weight is not applied.
777 </p></dd><dt><span class="term">
778 G4PSCellCharge
779 </span></dt><dd><p>
780 This class scored the total charge of particles which has stoped
781 in the volume.
782 </p></dd></dl></div><p>
783</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.Hits.G4VSDFil"></a>4.4.7. 
784<span class="emphasis"><em>G4VSDFilter</em></span> and its derived classes
785</h3></div></div></div><p>
786<span class="emphasis"><em>G4VSDFilter</em></span> is an abstract class that represents a track
787filter to be associated with <span class="emphasis"><em>G4VSensitiveDetector</em></span> or
788<span class="emphasis"><em>G4VPrimitiveScorer</em></span>. It defines a virtual method
789
790</p><div class="informalexample"><pre class="programlisting">
791 <span class="emphasis"><em>G4bool Accept(const G4Step*)</em></span>
792</pre></div><p>
793
794that should return <span class="emphasis"><em>true</em></span> if this particular step should be
795scored by the <span class="emphasis"><em>G4VSensitiveDetector</em></span> or
796<span class="emphasis"><em>G4VPrimitiveScorer</em></span>.
797</p><p>
798While the user can implement his/her own filter class, Geant4
799version 8.0 provides the following concrete filter classes:
800
801</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
802 G4SDChargedFilter
803 </span></dt><dd><p>
804 All charged particles are accepted.
805 </p></dd><dt><span class="term">
806 G4SDNeutralFilter
807 </span></dt><dd><p>
808 All neutral particles are accepted.
809 </p></dd><dt><span class="term">
810 G4SDParticleFilter
811 </span></dt><dd><p>
812 Particle species which are registered to this filter object by
813 <span class="emphasis"><em>Add("particle_name")</em></span> are accepted. More than one
814 species can be registered.
815 </p></dd><dt><span class="term">
816 G4SDKineticEnergyFilter
817 </span></dt><dd><p>
818 A track with kinetic energy greater than or equal to EKmin and
819 smaller than EKmin is accepted. EKmin and EKmax should be defined
820 as arguments of the constructor. The default values of EKmin and
821 EKmax are zero and DBL_MAX.
822 </p></dd><dt><span class="term">
823 G4SDParticleWithEnergyFilter
824 </span></dt><dd><p>
825 Combination of <span class="emphasis"><em>G4SDParticleFilter</em></span> and
826 <span class="emphasis"><em>G4SDParticleWithEnergyFilter</em></span>.
827 </p></dd></dl></div><p>
828</p><p>
829The use of the <span class="emphasis"><em>G4SDParticleFilter</em></span> class is demonstrated
830in <a href="ch04s04.html#programlist_Hits_2" title="Example 4.15. 
831An example of defining primitive sensitivity classes taken from
832ExN07DetectorConstruction.
833">Example 4.15</a>, where filters which accept gamma,
834electron, positron and electron/positron are defined.
835</p></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EvtBias.ScorImpRoul.Scor"></a>4.4.8. 
836Scoring for Event Biasing
837</h3></div></div></div><p>
838Scoring for Event Biasing (described in <a href="ch03s07.html" title="3.7. 
839Event Biasing Techniques
840">Section 3.7</a>) is a
841very specific use case whereby particle weights and fluxes through importance
842cells are required. The goals of the scoring technique are to:
843
844</p><div class="itemizedlist"><ul type="disc" compact><li><p>
845 appraise particle quantities related to special regions or
846 surfaces,
847 </p></li><li><p>
848 be applicable to all "cells" (physical volumes or replicas) of
849 a given geometry,
850 </p></li><li><p>
851 be customizable.
852 </p></li></ul></div><p>
853</p><p>
854Standard scoring must be provided for quantities such as tracks
855entering a cell, average weight of entering tracks, energy of
856entering tracks, and collisions inside the cell.
857</p><p>
858A number of scorers have been created for this specific appliction:
859</p><p>
860</p><div class="variablelist"><p class="title"><b></b></p><dl><dt><span class="term">
861 G4PSNofCollision
862 </span></dt><dd><p>
863 This scorer records the number of collisions that occur within a scored volume/cell.
864 There is the additional possibility to take into account the track weight
865 whilst scoring the number of collisions, via the following command:
866</p><div class="informalexample"><pre class="programlisting">
867 G4PSNofCollision* scorer1 = new G4PSNofCollision(psName="CollWeight");
868 scorer1-&gt;Weighted(true);
869</pre></div><p>
870 </p></dd><dt><span class="term">
871 G4PSPopulation
872 </span></dt><dd><p>
873 This scores the number of tracks within in a given cell per event.
874 </p></dd><dt><span class="term">
875 G4PSTrackLength
876 </span></dt><dd><p>
877 The track lengths within a cell are measured and if, additionally, the result is desired
878 to be weighted then the following code has to be implemented:
879</p><div class="informalexample"><pre class="programlisting">
880 G4PSTrackLength* scorer5 = new G4PSTrackLength(psName="SLW");
881 scorer5-&gt;Weighted(true);
882</pre></div><p>
883 Further if the energy track flux is required then the following should be
884 implemented:
885</p><div class="informalexample"><pre class="programlisting">
886 G4PSTrackLength* scorer6 = new G4PSTrackLength(psName="SLWE");
887 scorer6-&gt;Weighted(true);
888 scorer6-&gt;MultiplyKineticEnergy(true);
889 MFDet-&gt;RegisterPrimitive(scorer6);
890</pre></div><p>
891
892 Alternatively to measure the flux per unit velocity then:
893</p><div class="informalexample"><pre class="programlisting">
894 G4PSTrackLength* scorer7 = new G4PSTrackLength(psName="SLW_V");
895 scorer7-&gt;Weighted(true);
896 scorer7-&gt;DivideByVelocity(true);
897 MFDet-&gt;RegisterPrimitive(scorer7);
898</pre></div><p>
899
900 Finally to measure the flux energy per unit velocity then:
901</p><div class="informalexample"><pre class="programlisting">
902 G4PSTrackLength* scorer8 = new G4PSTrackLength(psName="SLWE_V");
903 scorer8-&gt;Weighted(true);
904 scorer8-&gt;MultiplyKineticEnergy(true);
905 scorer8-&gt;DivideByVelocity(true);
906 MFDet-&gt;RegisterPrimitive(scorer8);
907</pre></div><p>
908 </p></dd></dl></div><p>
909</p></div></div><div class="navfooter"><hr><table width="100%" summary="Navigation footer"><tr><td width="40%" align="left"><a accesskey="p" href="ch04s03.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><td width="20%" align="center"><a accesskey="u" href="ch04.html"><img src="AllResources/IconsGIF/up.gif" alt="Up"></a></td><td width="40%" align="right"> <a accesskey="n" href="ch04s05.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr><tr><td width="40%" align="left" valign="top">4.3. 
910Electromagnetic Field
911 </td><td width="20%" align="center"><a accesskey="h" href="index.html"><img src="AllResources/IconsGIF/home.gif" alt="Home"></a></td><td width="40%" align="right" valign="top"> 4.5. 
912Digitization
913</td></tr></table></div></body></html>
Note: See TracBrowser for help on using the repository browser.