source: trunk/documents/UserDoc/UsersGuides/ForApplicationDeveloper/html/Detector/hit.html

Last change on this file was 1208, checked in by garnier, 15 years ago

CVS update

File size: 29.5 KB
Line 
1<HTML>
2<TITLE>
3</TITLE>
4<!-- Changed by: Katsuya Amako, 21-Sep-1998 -->
5<!-- Proof read by: Joe Chuma,  30-Jun-1999 -->
6<!-- Changed by: Dennis Wright, 29-Nov-2001 -->
7<!-- Changed by: Makoto Asai,   30-Nov-2005 -->
8
9<BODY>
10<TABLE WIDTH="100%"><TR>
11<TD>
12
13
14<A HREF="index.html">
15<IMG SRC="../../../../resources/html/IconsGIF/Contents.gif" ALT="Contents"></A>
16<A HREF="electroMagneticField.html">
17<IMG SRC="../../../../resources/html/IconsGIF/Previous.gif" ALT="Previous"></A>
18<A HREF="digitization.html">
19<IMG SRC="../../../../resources/html/IconsGIF/Next.gif" ALT="Next"></A>
20</TD>
21<TD ALIGN="Right">
22<FONT SIZE="-1" COLOR="#238E23">
23<B>Geant4 User's Guide</B>
24<BR>
25<B>For Application Developers</B>
26<BR>
27<B>Detector Definition and Response</B>
28</FONT>
29</TD>
30</TR></TABLE>
31<BR><BR>
32
33<P ALIGN="Center">
34<FONT SIZE="+3" COLOR="#238E23">
35<B>4.4 Hits</B>
36</FONT>
37<BR><BR>
38
39<HR ALIGN="Center" SIZE="7%">
40<p>
41
42<a name="4.4.1">
43<H2>4.4.1 Hit</H2></a>
44
45A hit is a snapshot of the physical interaction of a track in the sensitive
46region of a detector. In it you can store information associated with a
47<i>G4Step</i> object. This information can be
48 <ul>
49  <li>the position and time of the step,
50  <li>the momentum and energy of the track,
51  <li>the energy deposition of the step,
52  <li>geometrical information,
53 </ul>
54or any combination of the above.
55<p>
56
57<b><i>G4VHit</i></b>
58<p>
59<i>G4VHit</i> is an abstract base class which represents a hit.  You must
60inherit this base class and derive your own concrete hit class(es).  The
61member data of your concrete hit class can be, and should be, your choice.
62</p>
63<p>
64<i>G4VHit</i> has two virtual methods, <tt>Draw()</tt> and <tt>Print()</tt>.
65To draw or print out your concrete hits, these methods should be implemented.
66How to define the drawing method is described in
67<a href="../Visualization/markertext.html">Section 8.9</a>.
68</p>
69
70<br>
71<b><i>G4THitsCollection</i></b>
72<p>
73<i>G4VHit</i> is an abstract class from which you derive your own concrete
74classes.  During the processing of a given event, represented by a
75<i>G4Event</i> object, many objects of the hit class will be produced,
76collected and associated with the event.  Therefore, for each concrete hit
77class you must also prepare a concrete class derived from
78<i>G4VHitsCollection</i>, an abstract class which represents a vector
79collection of user defined hits.
80</p>
81<p>
82<i>G4THitsCollection</i> is a template class derived from
83<i>G4VHitsCollection</i>, and the concrete hit collection class of a
84particular <i>G4VHit</i> concrete class can be instantiated from this template
85class.  Each object of a hit collection must have a unique name for each
86event.
87</p>
88<p>
89<i>G4Event</i> has a <i>G4HCofThisEvent</i> class object, that is a container
90class of collections of hits.  Hit collections are stored by their pointers,
91whose type is that of the base class.
92<p>
93
94<b>An example of a concrete hit class</b>
95<p>
96 Source listing 4.4.1 shows an example of a concrete hit class.
97<p>
98<center>
99<table border=2 cellpadding=10>
100<tr>
101<td>
102 <pre>
103#ifndef ExN04TrackerHit_h
104#define ExN04TrackerHit_h 1
105
106#include "G4VHit.hh"
107#include "G4THitsCollection.hh"
108#include "G4Allocator.hh"
109#include "G4ThreeVector.hh"
110
111class ExN04TrackerHit : public G4VHit
112{
113  public:
114
115      ExN04TrackerHit();
116      ~ExN04TrackerHit();
117      ExN04TrackerHit(const ExN04TrackerHit &right);
118      const ExN04TrackerHit& operator=(const ExN04TrackerHit &right);
119      int operator==(const ExN04TrackerHit &right) const;
120
121      inline void * operator new(size_t);
122      inline void operator delete(void *aHit);
123
124      void Draw() const;
125      void Print() const;
126
127  private:
128      G4double edep;
129      G4ThreeVector pos;
130
131  public:
132      inline void SetEdep(G4double de)
133      { edep = de; }
134      inline G4double GetEdep() const
135      { return edep; }
136      inline void SetPos(G4ThreeVector xyz)
137      { pos = xyz; }
138      inline G4ThreeVector GetPos() const
139      { return pos; }
140
141};
142
143typedef G4THitsCollection&lt;ExN04TrackerHit&gt; ExN04TrackerHitsCollection;
144
145extern G4Allocator&lt;ExN04TrackerHit&gt; ExN04TrackerHitAllocator;
146
147inline void* ExN04TrackerHit::operator new(size_t)
148{
149  void *aHit;
150  aHit = (void *) ExN04TrackerHitAllocator.MallocSingle();
151  return aHit;
152}
153
154inline void ExN04TrackerHit::operator delete(void *aHit)
155{
156  ExN04TrackerHitAllocator.FreeSingle((ExN04TrackerHit*) aHit);
157}
158
159#endif
160
161</pre>
162</td>
163</tr>
164<tr>
165<td align=center>
166 Source listing 4.4.1<BR>
167 An example of a concrete hit class.
168</td>
169</tr>
170</table></center>
171<p>
172<i>G4Allocator</i> is a class for fast allocation of objects to the heap
173through the paging mechanism.  For details of <i>G4Allocator</i>,
174refer to <a href="../Fundamentals/global.html">3.2.4</a>.  Use of
175<i>G4Allocator</i> is not mandatory, but it is recommended, especially for
176users who are not familiar with the C++ memory allocation
177mechanism or alternative tools of memory allocation.  On the other hand,
178note that <i>G4Allocator</i> is to be used <b>only</b> for the concrete class
179that is <b>not</b> used as a base class of any other classes.  For example,
180do <b>not</b> use the <i>G4Trajectory</i> class as a base class for a
181customized trajectory class, since <i>G4Trajectory</i> uses <i>G4Allocator</i>.
182</p>
183<p>
184
185<b>G4THitsMap</b>
186<p>
187<i>G4THitsMap</i> is an alternative to <i>G4THitsCollection</i>.
188<i>G4THitsMap</i> does not demand <i>G4VHit</i>, but instead any variable
189which can be mapped with an integer key.  Typically the key is a copy number
190of the volume, and the mapped value could for example be a double, such as
191the energy deposition in a volume. <i>G4THitsMap</i> is convenient for
192applications which do not need to output event-by-event data but instead just
193accumulate them.  All the <i>G4VPrimitiveScorer</i> classes discussed in
194<a href="#4.4.5">section 4.4.5</a> use <i>G4THitsMap</i>.
195</p>
196<p>
197<i>G4THitsMap</i> is derived from the <i>G4VHitsCollection</i> abstract base
198class and all objects of this class are also stored in <i>G4HCofThisEvent</i>
199at the end of an event. How to access a <i>G4THitsMap</i> object
200is discussed in the <a href="#4.4.5">following section</a>.
201
202<hr>
203<a name="4.4.2">
204<H2>4.4.2 Sensitive detector</H2></a>
205
206<b><i>G4VSensitiveDetector</i></b>
207<p>
208<i>G4VSensitiveDetector</i> is an abstract base class which represents
209a detector.  The principal mandate of a sensitive detector is the
210construction of hit objects using information from steps along a particle
211track.  The <tt>ProcessHits()</tt> method of <i>G4VSensitiveDetector</i>
212performs this task using <i>G4Step</i> objects as input.  In the case of a
213"Readout" geometry (see <a href="#4.4.3">Section 4.4.3</a>), objects of
214the <i>G4TouchableHistory</i> class may be used as an optional input.
215
216<p>
217Your concrete detector class should be instantiated with the unique name of
218your detector. The name can be associated with one or more global
219names with "/" as a delimiter for categorizing your detectors. For example
220 <pre>
221     myEMcal = new MyEMcal("/myDet/myCal/myEMcal");
222 </pre>
223where <tt>myEMcal</tt> is the name of your detector. The pointer to your
224sensitive detector must be set to one or more <i>G4LogicalVolume</i>
225objects to set the sensitivity of these volumes. The pointer should
226also be registered to <i>G4SDManager</i>, as described in
227<a href="#4.4.4">Section 4.4.4</a>.
228<p>
229 <i>G4VSensitiveDetector</i> has three major virtual methods.
230<p>
231 <dl>
232 <dt><tt>ProcessHits()</tt>
233  <dd>This method is invoked by <i>G4SteppingManager</i> when a step is
234  composed in the <i>G4LogicalVolume</i> which has the pointer to this
235  sensitive detector. The first argument of this method is a <i>G4Step</i>
236  object of the current step. The second argument is a
237  <i>G4TouchableHistory</i> object for the ``Readout geometry'' described
238  in the next section. The second argument is <tt>NULL</tt> if
239  ``Readout geometry'' is not assigned to this sensitive detector.
240  In this method, one or more <i>G4VHit</i> objects should be
241  constructed if the current step is meaningful for your detector.<p>
242 <dt><tt>Initialize()</tt>
243  <dd>This method is invoked at the beginning of each event. The argument
244  of this method is an object of the <i>G4HCofThisEvent</i> class. Hit
245  collections, where hits produced in this particular event are stored,
246  can be associated with the <i>G4HCofThisEvent</i> object in this
247  method. The hit collections associated with the <i>G4HCofThisEvent</i>
248  object during this method can be used for ``during the event processing''
249  digitization.<p>
250 <dt><tt>EndOfEvent()</tt>
251  <dd>This method is invoked at the end of each event. The argument of
252  this method is the same object as the previous method. Hit collections
253  occasionally created in your sensitive detector can be associated with
254  the <i>G4HCofThisEvent</i> object.
255 </dl>
256<p>
257
258<HR>
259<a name="4.4.3">
260<H2>4.4.3 Readout geometry</H2></a>
261
262This section describes how a ``Readout geometry'' can be defined.
263A Readout geometry is a virtual, parallel geometry for obtaining the
264channel number.
265<p>
266 As an example, the accordion calorimeter of <b>ATLAS</b> has a
267 complicated tracking geometry, however the readout can be done by
268 simple cylindrical sectors divided by theta, phi, and depth. Tracks
269 will be traced in the tracking geometry, the ``real'' one, and the
270 sensitive detector will have its own readout geometry Geant4 will
271 message to find to which ``readout'' cell the current hit belongs.
272 <center><img src="hit.src/RO.gif"></center>
273<p>
274The picture shows how this association is done in Geant4.  The first step is
275to associate a sensitive detector to a volume of the tracking geometry, in
276the usual way (see <a href="#4.4.2">Section 4.4.2</a>). The next step is to
277associate your <i>G4VReadoutGeometry</i> object to the sensitive detector.
278</p>
279<p>
280At tracking time, the base class <i>G4VReadoutGeometry</i> will provide to
281your sensitive detector code the <i>G4TouchableHistory</i> in the Readout
282geometry at the beginning of the step position (position of
283<i>PreStepPoint</i> of <i>G4Step</i>) and at this position only.
284</p>
285<p>
286This <i>G4TouchableHistory</i> is given to your sensitive detector code
287through the <i>G4VSensitiveDetector</i> virtual method:
288 <pre>
289        G4bool processHits(G4Step* aStep, G4TouchableHistory* ROhist);
290 </pre>
291 by the <tt>ROhist</tt> argument.
292<p>
293Thus, you will be able to use information from both the <i>G4Step</i> and
294the <i>G4TouchableHistory</i> coming from your Readout geometry.  Note that
295since the association is done through a sensitive detector object, it is
296perfectly possible to have several Readout geometries in parallel.
297<p>
298
299<b>Definition of a virtual geometry setup</b>
300<p>
301The base class for the implementation of a Readout geometry is
302<i>G4VReadoutGeometry</i>. This class has a single pure virtual protected
303method:
304 <pre>
305        virtual G4VPhysicalVolume* build() = 0;
306 </pre>
307which you must override in your concrete class. The <i>G4VPhysicalVolume</i> 
308pointer you will have to return is of the physical world of the Readout
309geometry.
310<p>
311The step by step procedure for constructing a Readout geometry is:
312 <ul>
313  <li>inherit from <i>G4VReadoutGeometry</i> to define a <i>MyROGeom</i> 
314      class;<p>
315  <li>implement the Readout geometry in the <tt>build()</tt> method,
316      returning the physical world of this geometry.<p>
317      The world is specified in the same way as for the
318      detector construction: a physical volume with no mother. The
319      axis system of this world is the same as the one of the
320      world for tracking.<p>
321      In this geometry you need to declare the sensitive parts
322      in the same way as in the tracking geometry: by setting a
323      non-<tt>NULL</tt> <i>G4VSensitiveDetector</i> pointer in, say, the
324      relevant <i>G4LogicalVolume</i> objects.  This sensitive class needs
325      to be there, but will not be used.<p>
326      You will also need to assign well defined materials to the volumes you
327      place in this geometry, but these materials are irrelevant since they
328      will not be seen by the tracking.  It is foreseen to allow the setting
329      of a <tt>NULL</tt> pointer in this case of the parallel geometry.<p>
330  <li>in the <tt>construct()</tt> method of your concrete
331      <i>G4VUserDetectorConstruction</i> class:<p>
332      <ul>
333        <li>instantiate your Readout geometry:
334            <pre>
335      MyROGeom* ROgeom = new MyROGeom("ROName");
336            </pre>
337        <li>build it:
338            <pre>
339      ROgeom-&gt;buildROGeometry();
340            </pre>
341            That will invoke your <tt>build()</tt> method.<p>
342        <li>Instantiate the sensitive detector which will
343            receive the <tt>ROGeom</tt> pointer, <tt>MySensitive</tt>, and
344            add this sensitive to the <i>G4SDManager</i>. Associate this
345            sensitive to the volume(s) of the tracking geometry as usual.<p>
346        <li>Associate the sensitive to the Readout geometry:
347            <pre>
348      MySensitive-&gt;SetROgeometry(ROgeom);
349            </pre>
350      </ul>
351 </ul>
352<p>
353
354<HR>
355<a name="4.4.4">
356<H2>4.4.4 <i>G4SDManager</i></H2></a>
357
358 <i>G4SDManager</i> is the singleton manager class for sensitive detectors.
359<p>
360
361<B>Activation / inactivation of sensitive detectors</B>
362<P>
363 The user interface commands <tt>activate</tt> and <tt>inactivate</tt>
364 are available to control your sensitive detectors. For example:
365 <p>
366 <tt>/hits/activate detector_name</tt><br>
367 <tt>/hits/inactivate detector_name</tt>
368 <p>
369 where <tt>detector_name</tt> can be the detector name or the category
370 name. For example, if your EM calorimeter is named <tt>/myDet/myCal/myEMcal</tt>,
371 <pre>
372     /hits/inactivate myCal
373 </pre>
374 will inactivate all detectors belonging to the <tt>myCal</tt> category.
375<p>
376
377<B>Access to the hit collections</B>
378<p>
379 Hit collections are accessed for various cases.
380 <p>
381 <ul>
382  <li>Digitization
383  <li>Event filtering in <i>G4VUserStackingAction</i>
384  <li>``End of event'' simple analysis
385  <li>Drawing / printing hits
386 </ul>
387 <p>
388The following is an example of how to access the hit collection of a
389particular concrete type:
390 <pre>
391  G4SDManager* fSDM = G4SDManager::GetSDMpointer();
392  G4RunManager* fRM = G4RunManager::GetRunManager();
393  G4int collectionID = fSDM-&gt;GetCollectionID("collection_name");
394  const G4Event* currentEvent = fRM-&gt;GetCurrentEvent();
395  G4HCofThisEvent* HCofEvent = currentEvent-&gt;GetHCofThisEvent();
396  <i>MyHitsCollection</i>* myCollection = (<i>MyHitsCollection</i>*)(HC0fEvent-&gt;GetHC(collectionID));
397 </pre>
398<p>
399
400<hr>
401<a name="4.4.5">
402<h2>4.4.5 <i>G4MultiFunctionalDetector</i> and
403<i>G4VPrimitiveScorer</i></h2></a>
404
405<i>G4MultiFunctionalDetector</i> is a concrete class derived from
406<i>G4VSensitiveDetector</i>.  Instead of implementing a user-specific detector
407class, <i>G4MultiFunctionalDetector</i> allows the user to register
408<i>G4VPrimitiveScorer</i> classes to build up the sensitivity.
409<i>G4MultiFunctionalDetector</i> should be instantiated in the users detector
410construction with its unique name and should be assigned to one or more
411<i>G4LogicalVolume</i>s.
412<p>
413
414<i>G4VPrimitiveScorer</i> is an abstract base class representing a class to be
415registered to <i>G4MultiFunctionalDetector</i> that creates a
416<i>G4THitsMap</i> object of one physics quantity for an event. Geant4 provides
417many concrete primitive scorer classes listed in the
418<a href="#4.4.6">next section</a>, and the user can also implement his/her
419own primitive scorers. Each primitive scorer object must be instantiated with
420a name that must be unique among primitive scorers registered in a
421<i>G4MultiFunctionalDetector</i>. Please note that a primitive scorer object
422must <b>not</b> be shared by more than one <i>G4MultiFunctionalDetector</i>
423object.
424<p>
425
426As mentioned in <a href="#4.4.1>4.4.1</a>, each <i>G4VPrimitiveScorer</i>
427generates one <i>G4THitsMap</i> object per event.  The name of the map object
428is the same as the name of the primitive scorer. Each of the concrete
429primitive scorers listed in the <a href="#4.4.6">next section</a>
430generates a <i>G4THitsMap&lt;G4double&gt;</i> that maps a <i>G4double</i> 
431value to its key integer number.  By default, the key is taken as the copy
432number of the <i>G4LogicalVolume</i> to which <i>G4MultiFunctionalDetector</i> 
433is assigned. In case the logical volume is uniquely placed in its mother
434volume and the mother is replicated, the copy number of its mother volume can
435be taken by setting the second argument of the <i>G4VPrimitiveScorer</i> 
436constructor, "<i>depth</i>" to 1, i.e. one level up. Furthermore, in case the
437key must consider more than one copy number of a different geometry hierarchy,
438the user can derive his/her own primitive scorer from the provided concrete
439class and implement the <i>GetIndex(G4Step*)</i> virtual method to return the
440unique key.
441
442<P>
443Source listing 4.4.2 shows an example of primitive sensitivity class
444definitions.
445<p>
446<center>
447<table border=2 cellpadding=10>
448<tr>
449<td>
450 <pre>
451void ExN07DetectorConstruction::SetupDetectors()
452{
453  G4String filterName, particleName;
454 
455  G4SDParticleFilter* gammaFilter = new G4SDParticleFilter(filterName="gammaFilter",particleName="gamma");
456  G4SDParticleFilter* electronFilter = new G4SDParticleFilter(filterName="electronFilter",particleName="e-");
457  G4SDParticleFilter* positronFilter = new G4SDParticleFilter(filterName="positronFilter",particleName="e+");
458  G4SDParticleFilter* epFilter = new G4SDParticleFilter(filterName="epFilter");
459  epFilter-&gt;add(particleName="e-");
460  epFilter-&gt;add(particleName="e+");
461   
462   
463  for(G4int i=0;i&lt;3;i++)
464  {
465   for(G4int j=0;j&lt;2;j++)
466   {
467    // Loop counter j = 0 : absorber
468    //                = 1 : gap
469    G4String detName = calName[i];
470    if(j==0)
471    { detName += "_abs"; }
472    else
473    { detName += "_gap"; }
474    G4MultiFunctionalDetector* det = new G4MultiFunctionalDetector(detName);
475 
476    //  The second argument in each primitive means the "level" of geometrical hierarchy,
477    // the copy number of that level is used as the key of the G4THitsMap.
478    //  For absorber (j = 0), the copy number of its own physical volume is used.
479    //  For gap (j = 1), the copy number of its mother physical volume is used, since there
480    // is only one physical volume of gap is placed with respect to its mother.
481    G4VPrimitiveScorer* primitive;
482    primitive = new G4PSEnergyDeposit("eDep",j);
483    det-&gt;RegisterPrimitive(primitive);
484    primitive = new G4PSNofSecondary("nGamma",j);
485    primitive-&gt;SetFilter(gammaFilter);
486    det-&gt;RegisterPrimitive(primitive);
487    primitive = new G4PSNofSecondary("nElectron",j);
488    primitive-&gt;SetFilter(electronFilter);
489    det-&gt;RegisterPrimitive(primitive);
490    primitive = new G4PSNofSecondary("nPositron",j);
491    primitive-&gt;SetFilter(positronFilter);
492    det-&gt;RegisterPrimitive(primitive);
493    primitive = new G4PSMinKinEAtGeneration("minEkinGamma",j);
494    primitive-&gt;SetFilter(gammaFilter);
495    det-&gt;RegisterPrimitive(primitive);
496    primitive = new G4PSMinKinEAtGeneration("minEkinElectron",j);
497    primitive-&gt;SetFilter(electronFilter);
498    det-&gt;RegisterPrimitive(primitive);
499    primitive = new G4PSMinKinEAtGeneration("minEkinPositron",j);
500    primitive-&gt;SetFilter(positronFilter);
501    det-&gt;RegisterPrimitive(primitive);
502    primitive = new G4PSTrackLength("trackLength",j);
503    primitive-&gt;SetFilter(epFilter);
504    det-&gt;RegisterPrimitive(primitive);
505    primitive = new G4PSNofStep("nStep",j);
506    primitive-&gt;SetFilter(epFilter);
507    det-&gt;RegisterPrimitive(primitive);
508 
509    G4SDManager::GetSDMpointer()-&gt;AddNewDetector(det);
510    if(j==0)
511    { layerLogical[i]-&gt;SetSensitiveDetector(det); }
512    else
513    { gapLogical[i]-&gt;SetSensitiveDetector(det); }
514   }
515  }
516}
517</pre>
518</td>
519</tr>
520<tr>
521<td align=center>
522 Source listing 4.4.2<BR>
523 An example of defining primitive sensitivity classes taken from <i>ExN07DetectorConstruction</i>.
524</td>
525</tr>
526</table></center>
527<p>
528
529Each <i>G4THitsMap</i> object can be accessed from <i>G4HCofThisEvent</i> with
530a unique collection ID number. This ID number can be obtained from
531<i>G4SDManager::GetCollectionID()</i> with a name of
532<i>G4MultiFunctionalDetector</i> and <i>G4VPrimitiveScorer</i> connected
533with a slush ("/").
534<i>G4THitsMap</i> has a [] operator taking the key value as an argument and
535returning <b>the pointer</b> of the value.  Please note that the [] operator
536returns <b>the pointer</b> of the value.  If you get zero from the []
537operator, it does <b>not</b> mean the value is zero, but that the provided key
538does not exist. The value itself is accessible with an astarisk ("*"). It is
539advised to check the validity of the returned pointer before accessing the
540value.  <i>G4THitsMap</i> also has a += operator in order to accumulate
541event data into run data. Source listing 4.4.3 shows the use of
542<i>G4THitsMap</i>.
543<p>
544<center>
545<table border=2 cellpadding=10>
546<tr>
547<td>
548 <pre>
549#include "ExN07Run.hh"
550#include "G4Event.hh"
551#include "G4HCofThisEvent.hh"
552#include "G4SDManager.hh"
553
554ExN07Run::ExN07Run()
555{
556  G4String detName[6] = {"Calor-A_abs","Calor-A_gap","Calor-B_abs","Calor-B_gap","Calor-C_abs","Calor-C_gap"};
557  G4String primNameSum[6] = {"eDep","nGamma","nElectron","nPositron","trackLength","nStep"};
558  G4String primNameMin[3] = {"minEkinGamma","minEkinElectron","minEkinPositron"};
559
560  G4SDManager* SDMan = G4SDManager::GetSDMpointer();
561  G4String fullName;
562  for(size_t i=0;i&lt;6;i++)
563  {
564    for(size_t j=0;j&lt;6;j++)
565    {
566      fullName = detName[i]+"/"+primNameSum[j];
567      colIDSum[i][j] = SDMan-&gt;GetCollectionID(fullName);
568    }
569    for(size_t k=0;k&lt;3;k++)
570    {
571      fullName = detName[i]+"/"+primNameMin[k];
572      colIDMin[i][k] = SDMan-&gt;GetCollectionID(fullName);
573    }
574  }
575}
576
577void ExN07Run::RecordEvent(const G4Event* evt)
578{
579  G4HCofThisEvent* HCE = evt-&gt;GetHCofThisEvent();
580  if(!HCE) return;
581  numberOfEvent++;
582  for(size_t i=0;i&lt;6;i++)
583  {
584    for(size_t j=0;j&lt;6;j++)
585    {
586      G4THitsMap&lt;G4double&gt;* evtMap = (G4THitsMap&lt;G4double&gt;*)(HCE-&gt;GetHC(colIDSum[i][j]));
587      mapSum[i][j] += *evtMap;
588    }
589    for(size_t k=0;k&lt;3;k++)
590    {
591      G4THitsMap&lt;G4double&gt;* evtMap = (G4THitsMap&lt;G4double&gt;*)(HCE-&gt;GetHC(colIDMin[i][k]));
592      std::map&lt;G4int,G4double*&gt;::iterator itr = evtMap-&gt;GetMap()-&gt;begin();
593      for(; itr != evtMap-&gt;GetMap()-&gt;end(); itr++)
594      {
595        G4int key = (itr-&gt;first);
596        G4double val = *(itr-&gt;second);
597        G4double* mapP = mapMin[i][k][key];
598        if( mapP && (val&gt;*mapP) ) continue;
599        mapMin[i][k].set(key,val);
600      }
601    }
602  }
603}
604</pre>
605</td>
606</tr>
607<tr>
608<td align=center>
609 Source listing 4.4.3<BR>
610 An example of accessing to <i>G4THitsMap</i> objects.
611</td>
612</tr>
613</table></center>
614<p>
615 
616<HR>
617<a name="4.4.6">
618<H2>4.4.6 Concrete classes of <i>G4VPrimitiveScorer</i></H2></a>
619
620With Geant4 version 8.0, several concrete primitive scorer classes are
621provided, all of which are derived from the <i>G4VPrimitiveScorer</i> abstract
622base class and which are to be registered to <i>G4MultiFunctionalDetector</i>.
623Each of them contains one <i>G4THitsMap</i> object and scores a simple double
624value for each key.
625
626<p>
627<b>Track length scorers</b><p>
628
629<DL>
630<DT>G4PSTrackLength</DT>
631<DD>The track length is defined as the sum of step lengths of the particles
632inside the cell. Bt default, the track weight is not taken into account,
633but could be used as a multiplier of each step length if the
634<i>Weighted()</i> method of this class object is invoked.</DD>
635<p>
636<DT>G4PSPassageTrackLength</DT>
637<DD>The passage track length is the same as the track length in
638<i>G4PSTrackLength</i>, except that only tracks which pass through the volume
639are taken into account. It means newly-generated or stopped tracks inside the
640cell are excluded from the calculation. By default, the track weight is not 
641taken into account, but could be used as a multiplier of each step length if
642the <i>Weighted()</i> method of this class object is invoked.</DD>
643</DL>
644
645<p>
646<b>Deposited energy scorers</b><p>
647
648<DL>
649<DT>G4PSEnergyDeposit</DT>
650<DD>This scorer stores a sum of particles' energy deposits at each step in the
651cell. The particle weight is multiplied at each step.</DD>
652<p>
653<DT>G4PSDoseDeposit</DT>
654<DD>In some cases, dose is a more convenient way to evaluate the effect of
655energy deposit in a cell than simple deposited energy.  The dose deposit is
656defined by the sum of energy deposits at each step in a cell divided by the
657mass of the cell. The mass is calculated from the density and volume of the
658cell taken from the methods of <i>G4VSolid</i> and <i>G4LogicalVolume</i>.
659The particle weight is multiplied at each step.</DD>
660</DL>
661
662<p>
663<b>Current and flux scorers</b><p>
664There are two different definitions of a particle's flow for a given geometry.
665One is a current and the other is a flux.  In our scorers, the current is
666simply defined as the number of particles (with the particle's weight) at a
667certain surface or volume, while the flux takes the particle's injection angle
668to the geometry into account.  The current and flux are usually defined at a
669surface, but volume current and volume flux are also provided.
670<p>
671<DL>
672<DT>G4PSFlatSurfaceCurrent</DT>
673<DD>Flat surface current is a surface based scorer. The present implementation
674is limited to scoring only at the -Z surface of a <i>G4Box</i> solid. The
675quantity is defined by the number of tracks that reach the surface.
676The user must choose a direction of the particle to be scored.  The choices
677are fCurrent_In, fCurrent_Out, or fCurrent_InOut, one of which must be entered
678as the second argument of the constructor. Here, fCurrent_In scores incoming
679particles to the cell, while fCurrent_Out scores only outgoing particles from
680the cell.  fCurrent_InOut scores both directions. The current is multiplied by
681particle weight and is normalized for a unit area.</DD><p>
682<DT>G4PSSphereSurfaceCurrent</DT>
683<DD>Sphere surface current is a surface based scorer, and similar to the
684G4PSFlatSurfaceCurrent. The only difference is that the surface is defined at
685the <b>inner surface</b> of a G4Sphere solid.</DD><p>
686<DT>G4PSPassageCurrent</DT>
687<DD>Passage current is a volume-based scorer. The current is defined by the
688number of tracks that pass through the volume. A particle weight is applied at
689the exit point.  A passage current is defined for a volume.</DD><p>
690
691<DT>G4PSFlatSurfaceFlux</DT>
692<DD>Flat surface flux is a surface based flux scorer. The surface flux is
693defined by the number of tracks that reach the surface. The expression of
694surface flux is given by the sum of W/cos(t)/A, where W, t and A represent
695particle weight, injection angle of particle with respect to the surface
696normal, and area of the surface. The user must enter one of the particle
697directions, fFlux_In, fFlux_Out, or fFlux_InOut in the constructor.
698Here, fFlux_In scores incoming particles to the cell, while fFlux_Out scores
699outgoing particles from the cell.  fFlux_InOut scores both directions.</DD><p>
700<DT>G4PSCellFlux</DT>
701<DD>Cell flux is a volume based flux scorer. The cell flux is defined by a
702track length (L) of the particle inside a volume divided by the volume (V) of
703this cell. The track length is calculated by a sum of the step lengths in the
704cell. The expression for cell flux is given by the sum of (W*L)/V, where W is
705a particle weight, and is multiplied by the track length at each step.</DD><p>
706<DT>G4PSPassageCellFlux</DT>
707<DD>Passage cell flux is a volume based scorer similar to <i>G4PSCellFlux</i>.
708The only difference is that tracks which pass through a cell are taken into
709account. It means generated or stopped tracks inside the volume are excluded
710from the calculation.</DD>
711</DL>
712
713<p>
714<b>Other scorers</b><p>
715<DL>
716<DT>G4PSMinKinEAtGeneration</DT>
717<DD>This scorer records the minimum kinetic energy of secondary particles at
718their production point in the volume in an event. This primitive scorer does
719not integrate the quantity, but records the minimum quantity. </DD><p>
720<DT>G4PSNofSecondary</DT>
721<DD>This class scores the number of secondary particles generated in the
722volume. The weight of the secondary track is taken into account. </DD><p>
723<DT>G4PSNofStep</DT>
724<DD>This class scores the number of steps in the cell. A particle weight is
725    not applied.</DD>
726<DT>G4PSCellCharge</DT>
727<DD>This class scored the total charge of particles which has stoped in the
728volume.</DD>
729</DL>
730
731<hr>
732<a name="4.4.7">
733<H2>4.4.7 <i>G4VSDFilter</i> and its derived classes</H2></a>
734
735<i>G4VSDFilter</i> is an abstract class that represents a track filter to be
736associated with <i>G4VSensitiveDetector</i> or <i>G4VPrimitiveScorer</i>. It
737defines a virtual method
738<DL><DD><i>G4bool Accept(const G4Step*)</i></DD></DL>
739that should return <i>true</i> if this particular step should be scored by
740the <i>G4VSensitiveDetector</i> or <i>G4VPrimitiveScorer</i>.
741
742<p>
743While the user can implement his/her own filter class, Geant4 version 8.0
744provides the following concrete filter classes:
745<p>
746<DL>
747<DT>G4SDChargedFilter</DL>
748<DD>All charged particles are accepted.</DD><p>
749<DT>G4SDNeutralFilter</DL>
750<DD>All neutral particles are accepted.</DD><p>
751<DT>G4SDParticleFilter</DL>
752<DD>Particle species which are registered to this filter object by
753   <i>Add("particle_name")</i> are accepted.  More than one species can be
754   registered.</DD><p>
755<DT>G4SDKineticEnergyFilter</DT>
756<DD>A track with kinetic energy greater than or equal to EKmin and smaller
757    than EKmin is accepted. EKmin and EKmax should be defined as arguments
758    of the constructor.  The default values of EKmin and EKmax are zero and
759    DBL_MAX.</DD><p>
760<DT>G4SDParticleWithEnergyFilter</DT>
761<DD>Combination of <i>G4SDParticleFilter</i> and
762   <i>G4SDParticleWithEnergyFilter</i>.</DD><p>
763</DL>
764
765The use of the <i>G4SDParticleFilter</i> class is demonstrated in the Source
766Listing 4.4.2, where filters which accept gamma, electron, positron and
767electron/positron are defined.
768
769<br><br>
770<hr>
771<a href="../../../../Authors/html/subjectsToAuthors.html">
772<i>About the authors</i></a>
773
774</body>
775</html>
Note: See TracBrowser for help on using the repository browser.