| 1 | <html><head><meta http-equiv="Content-Type" content="text/html; charset=ISO-8859-1"><title>3.7. Event Biasing Techniques</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="ch03.html" title="Chapter 3. Toolkit Fundamentals"><link rel="prev" href="ch03s06.html" title="3.6. Event Generator Interface"><link rel="next" href="ch04.html" title="Chapter 4. Detector Definition and Response"><script language="JavaScript">
|
|---|
| 2 | function 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">3.7.
|
|---|
| 9 | Event Biasing Techniques
|
|---|
| 10 | </th></tr><tr><td width="20%" align="left"><a accesskey="p" href="ch03s06.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><th width="60%" align="center">Chapter 3.
|
|---|
| 11 | Toolkit Fundamentals
|
|---|
| 12 | </th><td width="20%" align="right"> <a accesskey="n" href="ch04.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.EvtBias"></a>3.7.
|
|---|
| 13 | Event Biasing Techniques
|
|---|
| 14 | </h2></div></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EvtBias.ScorImpRoul"></a>3.7.1.
|
|---|
| 15 | Scoring, Geometrical Importance Sampling and Weight Roulette
|
|---|
| 16 | </h3></div></div></div><p>
|
|---|
| 17 | Geant4 provides event biasing techniques which may be used to save
|
|---|
| 18 | computing time in such applications as the simulation of radiation
|
|---|
| 19 | shielding. These are <span class="emphasis"><em>geometrical splitting
|
|---|
| 20 | </em></span> and <span class="emphasis"><em>Russian roulette</em></span>
|
|---|
| 21 | (also called geometrical importance sampling), and
|
|---|
| 22 | <span class="emphasis"><em>weight roulette</em></span>. Scoring is carried out by <span class="emphasis"><em>G4MultiFunctionalDetector</em></span> (see <a href="ch04s04.html#sect.Hits.G4Multi" title="4.4.5.
|
|---|
| 23 | G4MultiFunctionalDetector and
|
|---|
| 24 | G4VPrimitiveScorer
|
|---|
| 25 | ">Section 4.4.5</a> and <a href="ch04s04.html#sect.Hits.G4VPrim" title="4.4.6.
|
|---|
| 26 | Concrete classes of G4VPrimitiveScorer
|
|---|
| 27 | ">Section 4.4.6</a>) using the standard Geant4 scoring technique.
|
|---|
| 28 | Biasing specific scorers have been implemented and are described within
|
|---|
| 29 | <span class="emphasis"><em>G4MultiFunctionDetector</em></span> documentation. In this chapter, it is assumed that
|
|---|
| 30 | the reader is familiar with both the usage of Geant4 and the
|
|---|
| 31 | concepts of importance sampling. More detailed documentation may be
|
|---|
| 32 | found in the documents
|
|---|
| 33 | <a href="http://geant4.web.cern.ch/geant4/collaboration/working_groups/geometry/index.shtml" target="_top">
|
|---|
| 34 | 'Latest development in importance sampling and scoring'
|
|---|
| 35 | </a>.
|
|---|
| 36 | A detailed description of different use-cases which employ the sampling
|
|---|
| 37 | and scoring techniques can be found in the document
|
|---|
| 38 | <a href="http://geant4.web.cern.ch/geant4/collaboration/working_groups/geometry/index.shtml" target="_top">
|
|---|
| 39 | 'Scoring and geometrical importance sampling use cases'
|
|---|
| 40 | </a>.
|
|---|
| 41 | </p><p>
|
|---|
| 42 | The purpose of importance sampling is to save computing time by
|
|---|
| 43 | sampling less often the particle histories entering "less
|
|---|
| 44 | important" geometry regions, and more often in more "important"
|
|---|
| 45 | regions. Given the same amount of computing time, an
|
|---|
| 46 | importance-sampled and an analogue-sampled simulation must show
|
|---|
| 47 | equal mean values, while the importance-sampled simulation will
|
|---|
| 48 | have a decreased variance.
|
|---|
| 49 | </p><p>
|
|---|
| 50 | The implementation of scoring is independent of the implementation
|
|---|
| 51 | of importance sampling. However both share common concepts.
|
|---|
| 52 | <span class="emphasis"><em>Scoring and importance sampling apply to particle types chosen
|
|---|
| 53 | by the user</em></span>, which should be borne in mind when interpreting the
|
|---|
| 54 | output of any biased simulation.
|
|---|
| 55 | </p><p>
|
|---|
| 56 | Examples on how to use scoring and importance sampling may be found
|
|---|
| 57 | in <code class="literal">examples/extended/biasing</code> and
|
|---|
| 58 | <code class="literal">examples/advanced/Tiara</code>.
|
|---|
| 59 | </p><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.Geom"></a>3.7.1.1.
|
|---|
| 60 | Geometries
|
|---|
| 61 | </h4></div></div></div><p>
|
|---|
| 62 | The kind of scoring referred to in this note and the importance
|
|---|
| 63 | sampling apply to spatial cells provided by the user.
|
|---|
| 64 | </p><p>
|
|---|
| 65 | A <span class="bold"><strong>cell</strong></span> is a physical volume (further specified
|
|---|
| 66 | by it's replica number, if the volume is a replica). Cells may be defined
|
|---|
| 67 | in two kinds of geometries:
|
|---|
| 68 |
|
|---|
| 69 | </p><div class="orderedlist"><ol type="1" compact><li><p>
|
|---|
| 70 | <span class="bold"><strong>mass geometry</strong></span>: the geometry setup of the
|
|---|
| 71 | experiment to be simulated. Physics processes apply to this geometry.
|
|---|
| 72 | </p></li><li><p>
|
|---|
| 73 | <span class="bold"><strong>parallel-geometry</strong></span>: a geometry constructed
|
|---|
| 74 | to define the physical volumes according to which scoring and/or importance
|
|---|
| 75 | sampling is applied.
|
|---|
| 76 | </p></li></ol></div><p>
|
|---|
| 77 | </p><p>
|
|---|
| 78 | The user has the choice to score and/or sample by importance the
|
|---|
| 79 | particles of the chosen type, according to mass geometry or to
|
|---|
| 80 | parallel geometry. It is possible to utilize several parallel
|
|---|
| 81 | geometries in addition to the mass geometry. This provides the user
|
|---|
| 82 | with a lot of flexibility to define separate geometries for
|
|---|
| 83 | different particle types in order to apply scoring or/and
|
|---|
| 84 | importance sampling.
|
|---|
| 85 | </p><p>
|
|---|
| 86 | </p><div class="note" style="margin-left: 0.5in; margin-right: 0.5in;"><h3 class="title">Note</h3>
|
|---|
| 87 | Parallel geometries should be constructed using the implementation as
|
|---|
| 88 | described in <a href="ch04s07.html" title="4.7.
|
|---|
| 89 | Parallel Geometries
|
|---|
| 90 | ">Section 4.7</a>.
|
|---|
| 91 |
|
|---|
| 92 | There are a few conditions for parallel geometries:
|
|---|
| 93 |
|
|---|
| 94 | <div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 95 | The world volume for parallel and mass geometries must be identical copies.
|
|---|
| 96 | </p></li><li><p>
|
|---|
| 97 | Scoring and importance cells must not share boundaries with the world volume.
|
|---|
| 98 | </p></li></ul></div></div><p>
|
|---|
| 99 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.ChgSamp"></a>3.7.1.2.
|
|---|
| 100 | Changing the Sampling
|
|---|
| 101 | </h4></div></div></div><p>
|
|---|
| 102 | Samplers are higher level tools which perform the necessary
|
|---|
| 103 | changes of the Geant4 sampling in order to apply importance
|
|---|
| 104 | sampling and weight roulette.
|
|---|
| 105 | </p><p>
|
|---|
| 106 | Variance reduction (and scoring through the
|
|---|
| 107 | <span class="emphasis"><em>G4MultiFunctionalDetector</em></span>)
|
|---|
| 108 | may be combined arbitrarily for chosen particle types and may be applied to the
|
|---|
| 109 | mass or to parallel geometries.
|
|---|
| 110 | </p><p>
|
|---|
| 111 | The <code class="literal">G4GeometrySampler</code> can be applied equally to mass or
|
|---|
| 112 | parallel geometries with an abstract interface supplied by <code class="literal">G4VSampler</code>.
|
|---|
| 113 | <code class="literal">G4VSampler</code> provides
|
|---|
| 114 | <code class="literal">Prepare...</code> methods and a <code class="literal">Configure</code>
|
|---|
| 115 | method:
|
|---|
| 116 |
|
|---|
| 117 | <a name="anchor_EvtBias_G4VSampler"></a>
|
|---|
| 118 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 119 | class G4VSampler
|
|---|
| 120 | {
|
|---|
| 121 | public:
|
|---|
| 122 | G4VSampler();
|
|---|
| 123 | virtual ~G4VSampler();
|
|---|
| 124 | virtual void PrepareImportanceSampling(G4VIStore *istore,
|
|---|
| 125 | const G4VImportanceAlgorithm
|
|---|
| 126 | *ialg = 0) = 0;
|
|---|
| 127 | virtual void PrepareWeightRoulett(G4double wsurvive = 0.5,
|
|---|
| 128 | G4double wlimit = 0.25,
|
|---|
| 129 | G4double isource = 1) = 0;
|
|---|
| 130 | virtual void PrepareWeightWindow(G4VWeightWindowStore *wwstore,
|
|---|
| 131 | G4VWeightWindowAlgorithm *wwAlg = 0,
|
|---|
| 132 | G4PlaceOfAction placeOfAction =
|
|---|
| 133 | onBoundary) = 0;
|
|---|
| 134 | virtual void Configure() = 0;
|
|---|
| 135 | virtual void ClearSampling() = 0;
|
|---|
| 136 | virtual G4bool IsConfigured() const = 0;
|
|---|
| 137 | };
|
|---|
| 138 | </pre></div><p>
|
|---|
| 139 | </p><p>
|
|---|
| 140 | The methods for setting up the desired combination need specific
|
|---|
| 141 | information:
|
|---|
| 142 |
|
|---|
| 143 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 144 | Importance sampling: message <code class="literal">PrepareImportanceSampling</code>
|
|---|
| 145 | with a
|
|---|
| 146 | <a href="ch03s07.html#anchor_EvtBias_G4VIStore">
|
|---|
| 147 | <code class="literal">G4VIStore</code></a>
|
|---|
| 148 | and optionally a <code class="literal">G4VImportanceAlgorithm</code>
|
|---|
| 149 | </p></li><li><p>
|
|---|
| 150 | Weight window: message <code class="literal">PrepareWeightWindow</code> with the
|
|---|
| 151 | arguments:
|
|---|
| 152 | </p><div class="itemizedlist"><ul type="circle" compact><li><p>
|
|---|
| 153 | <span class="emphasis"><em>*wwstore</em></span>: a
|
|---|
| 154 | <code class="literal">G4VWeightWindowStore</code> for retrieving the lower
|
|---|
| 155 | weight bounds for the energy-space cells
|
|---|
| 156 | </p></li><li><p>
|
|---|
| 157 | <span class="emphasis"><em>*wwAlg</em></span>: a
|
|---|
| 158 | <code class="literal">G4VWeightWindowAlgorithm</code> if a customized algorithm
|
|---|
| 159 | should be used
|
|---|
| 160 | </p></li><li><p>
|
|---|
| 161 | <span class="emphasis"><em>placeOfAction</em></span>: a
|
|---|
| 162 | <code class="literal">G4PlaceOfAction</code> specifying where to perform the
|
|---|
| 163 | biasing
|
|---|
| 164 | </p></li></ul></div><p>
|
|---|
| 165 | </p></li><li><p>
|
|---|
| 166 | Weight roulette: message <code class="literal">PrepareWeightRoulett</code> with the
|
|---|
| 167 | optional parameters:
|
|---|
| 168 | </p><div class="itemizedlist"><ul type="circle" compact><li><p>
|
|---|
| 169 | <span class="emphasis"><em>wsurvive</em></span>: survival weight
|
|---|
| 170 | </p></li><li><p>
|
|---|
| 171 | <span class="emphasis"><em>wlimit</em></span>: minimal allowed
|
|---|
| 172 | value of weight * source importance / cell importance
|
|---|
| 173 | </p></li><li><p>
|
|---|
| 174 | <span class="emphasis"><em>isource</em></span>: importance of the source cell
|
|---|
| 175 | </p></li></ul></div><p>
|
|---|
| 176 | </p></li></ul></div><p>
|
|---|
| 177 | </p><p>
|
|---|
| 178 | Each object of a sampler class is responsible for one particle
|
|---|
| 179 | type. The particle type is given to the constructor of the sampler
|
|---|
| 180 | classes via the particle type name, e.g. "neutron". Depending on
|
|---|
| 181 | the specific purpose, the <code class="literal">Configure()</code> of a sampler will
|
|---|
| 182 | set up specialized processes (derived from <code class="literal">G4VProcess</code>) for
|
|---|
| 183 | transportation in the parallel geometry, importance
|
|---|
| 184 | sampling and weight roulette for the given particle type. When
|
|---|
| 185 | <code class="literal">Configure()</code> is invoked the sampler places the processes in
|
|---|
| 186 | the correct order independent of the order in which user invoked
|
|---|
| 187 | the <code class="literal">Prepare...</code> methods.
|
|---|
| 188 | </p><p>
|
|---|
| 189 | </p><div class="note" style="margin-left: 0.5in; margin-right: 0.5in;"><h3 class="title">Note</h3><p>
|
|---|
| 190 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 191 | The <code class="literal">Prepare...()</code> functions may each only be invoked
|
|---|
| 192 | once.
|
|---|
| 193 | </p></li><li><p>
|
|---|
| 194 | To configure the sampling the function <code class="literal">Configure()</code>
|
|---|
| 195 | must be called <span class="emphasis"><em>after</em></span> the
|
|---|
| 196 | <code class="literal">G4RunManager</code> has been initialized and the PhysicsList has
|
|---|
| 197 | been instantiated.
|
|---|
| 198 | </p></li></ul></div><p>
|
|---|
| 199 | </p></div><p>
|
|---|
| 200 | </p><p>
|
|---|
| 201 | The interface and framework are demonstrated in the
|
|---|
| 202 | <code class="literal">examples/extended/biasing</code> directory, with the main changes being to the
|
|---|
| 203 | G4GeometrySampler class and the fact that in the parallel case the WorldVolume
|
|---|
| 204 | is a copy of the Mass World.
|
|---|
| 205 |
|
|---|
| 206 | The parallel geometry now has to inherit from
|
|---|
| 207 | <span class="emphasis"><em>G4VUserParallelWorld</em></span> which also has the <span class="emphasis"><em>GetWorld()</em></span> method
|
|---|
| 208 | in order to retrieve a copy of the mass geometry WorldVolume.
|
|---|
| 209 |
|
|---|
| 210 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 211 | class B02ImportanceDetectorConstruction : public G4VUserParallelWorld
|
|---|
| 212 | ghostWorld = GetWorld();
|
|---|
| 213 | </pre></div><p>
|
|---|
| 214 |
|
|---|
| 215 | </p><p>
|
|---|
| 216 | The constructor for <span class="emphasis"><em>G4GeometrySampler</em></span> takes a pointer to
|
|---|
| 217 | the physical world volume and the particle type name (e.g. "neutron") as arguments.
|
|---|
| 218 | In a single mass geometry the sampler is created as follows:
|
|---|
| 219 |
|
|---|
| 220 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 221 | G4GeometrySampler mgs(detector->GetWorldVolume(),"neutron");
|
|---|
| 222 | mgs.SetParallel(false);
|
|---|
| 223 | </pre></div><p>
|
|---|
| 224 |
|
|---|
| 225 |
|
|---|
| 226 | Whilst the following lines of code are required in order to set up the sampler for the
|
|---|
| 227 | parallel geometry case:
|
|---|
| 228 |
|
|---|
| 229 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 230 | G4VPhysicalVolume* ghostWorld = pdet->GetWorldVolume();
|
|---|
| 231 |
|
|---|
| 232 | G4GeometrySampler pgs(ghostWorld,"neutron");
|
|---|
| 233 |
|
|---|
| 234 | pgs.SetParallel(true);
|
|---|
| 235 | </pre></div><p>
|
|---|
| 236 |
|
|---|
| 237 |
|
|---|
| 238 | Also note that the preparation and configuration of the samplers has to be
|
|---|
| 239 | carried out <span class="emphasis"><em>after</em></span> the instantiation of the UserPhysicsList
|
|---|
| 240 | and after the initialisation of the <span class="emphasis"><em>G4RunManager</em></span>:
|
|---|
| 241 |
|
|---|
| 242 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 243 | pgs.PrepareImportanceSampling(&aIstore, 0);
|
|---|
| 244 | pgs.Configure();
|
|---|
| 245 | </pre></div><p>
|
|---|
| 246 |
|
|---|
| 247 | Due to the fact that biasing is a process and has to be inserted after all the
|
|---|
| 248 | other processes have been created.
|
|---|
| 249 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.ImpSamp"></a>3.7.1.3.
|
|---|
| 250 | Importance Sampling
|
|---|
| 251 | </h4></div></div></div><p>
|
|---|
| 252 | Importance sampling acts on particles crossing boundaries
|
|---|
| 253 | between "importance cells". The action taken depends on the
|
|---|
| 254 | importance values assigned to the cells. In general a particle
|
|---|
| 255 | history is either split or Russian roulette is played if the
|
|---|
| 256 | importance increases or decreases, respectively. A weight assigned
|
|---|
| 257 | to the history is changed according to the action taken.
|
|---|
| 258 | </p><p>
|
|---|
| 259 | The tools provided for importance sampling require the user to
|
|---|
| 260 | have a good understanding of the physics in the problem. This is
|
|---|
| 261 | because the user has to decide which particle types require
|
|---|
| 262 | importance sampled, define the cells, and assign importance values
|
|---|
| 263 | to the cells. If this is not done properly the results cannot be
|
|---|
| 264 | expected to describe a real experiment.
|
|---|
| 265 | </p><p>
|
|---|
| 266 | The assignment of importance values to a cell is done using an
|
|---|
| 267 | importance store described below.
|
|---|
| 268 | </p><p>
|
|---|
| 269 | An "importance store" with the interface <code class="literal">G4VIStore</code> is
|
|---|
| 270 | used to store importance values related to cells. In order to do
|
|---|
| 271 | importance sampling the user has to create an object (e.g. of class
|
|---|
| 272 | <code class="literal">G4IStore</code>) of type <code class="literal">G4VIStore</code>. The samplers may be
|
|---|
| 273 | given a <code class="literal">G4VIStore</code>. The user fills the store with cells and
|
|---|
| 274 | their importance values.
|
|---|
| 275 | </p><p>
|
|---|
| 276 | An importance store has to be constructed with a reference to
|
|---|
| 277 | the world volume of the geometry used for importance sampling. This
|
|---|
| 278 | may be the world volume of the mass or of a parallel geometry.
|
|---|
| 279 | Importance stores derive from the interface <code class="literal">G4VIStore</code>:
|
|---|
| 280 |
|
|---|
| 281 | <a name="anchor_EvtBias_G4VIStore"></a>
|
|---|
| 282 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 283 | class G4VIStore
|
|---|
| 284 | {
|
|---|
| 285 | public:
|
|---|
| 286 | G4VIStore();
|
|---|
| 287 | virtual ~G4VIStore();
|
|---|
| 288 | virtual G4double GetImportance(const G4GeometryCell &gCell) const = 0;
|
|---|
| 289 | virtual G4bool IsKnown(const G4GeometryCell &gCell) const = 0;
|
|---|
| 290 | virtual const G4VPhysicalVolume &GetWorldVolume() const = 0;
|
|---|
| 291 | };
|
|---|
| 292 | </pre></div><p>
|
|---|
| 293 | </p><p>
|
|---|
| 294 | A concrete implementation of an importance store is provided by
|
|---|
| 295 | the class <code class="literal">G4VStore</code>. The <span class="emphasis"><em>public</em></span>
|
|---|
| 296 | part of the class is:
|
|---|
| 297 |
|
|---|
| 298 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 299 | class G4IStore : public G4VIStore
|
|---|
| 300 | {
|
|---|
| 301 | public:
|
|---|
| 302 | explicit G4IStore(const G4VPhysicalVolume &worldvolume);
|
|---|
| 303 | virtual ~G4IStore();
|
|---|
| 304 | virtual G4double GetImportance(const G4GeometryCell &gCell) const;
|
|---|
| 305 | virtual G4bool IsKnown(const G4GeometryCell &gCell) const;
|
|---|
| 306 | virtual const G4VPhysicalVolume &GetWorldVolume() const;
|
|---|
| 307 | void AddImportanceGeometryCell(G4double importance,
|
|---|
| 308 | const G4GeometryCell &gCell);
|
|---|
| 309 | void AddImportanceGeometryCell(G4double importance,
|
|---|
| 310 | const G4VPhysicalVolume &,
|
|---|
| 311 | G4int aRepNum = 0);
|
|---|
| 312 | void ChangeImportance(G4double importance,
|
|---|
| 313 | const G4GeometryCell &gCell);
|
|---|
| 314 | void ChangeImportance(G4double importance,
|
|---|
| 315 | const G4VPhysicalVolume &,
|
|---|
| 316 | G4int aRepNum = 0);
|
|---|
| 317 | G4double GetImportance(const G4VPhysicalVolume &,
|
|---|
| 318 | G4int aRepNum = 0) const ;
|
|---|
| 319 | private: .....
|
|---|
| 320 | };
|
|---|
| 321 | </pre></div><p>
|
|---|
| 322 | </p><p>
|
|---|
| 323 | The member function <code class="literal">AddImportanceGeometryCell()</code> enters
|
|---|
| 324 | a cell and an importance value into the importance store. The
|
|---|
| 325 | importance values may be returned either according to a physical
|
|---|
| 326 | volume and a replica number or according to a
|
|---|
| 327 | <code class="literal">G4GeometryCell</code>. The user must be aware of the
|
|---|
| 328 | interpretation of assigning importance values to a cell.
|
|---|
| 329 | If scoring is also implemented then this is attached to logical volumes, in
|
|---|
| 330 | which case the physical volume and replica number method should be used for
|
|---|
| 331 | assigning importance values. See <code class="literal">examples/extended/biasing
|
|---|
| 332 | B01</code> and <code class="literal">B02</code> for examples of this.
|
|---|
| 333 | </p><p>
|
|---|
| 334 | </p><div class="note" style="margin-left: 0.5in; margin-right: 0.5in;"><h3 class="title">Note</h3><p>
|
|---|
| 335 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 336 | An importance value must be assigned to every cell.
|
|---|
| 337 | </p></li></ul></div><p>
|
|---|
| 338 | </p></div><p>
|
|---|
| 339 | </p><p>
|
|---|
| 340 | The different cases:
|
|---|
| 341 |
|
|---|
| 342 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 343 | <span class="emphasis"><em>Cell is not in store</em></span>
|
|---|
| 344 | </p><p>
|
|---|
| 345 | Not filling a certain cell in the store will cause an
|
|---|
| 346 | exception.
|
|---|
| 347 | </p><p>
|
|---|
| 348 | </p></li><li><p>
|
|---|
| 349 | <span class="emphasis"><em>Importance value = zero</em></span>
|
|---|
| 350 | </p><p>
|
|---|
| 351 | Tracks of the chosen particle type will be killed.
|
|---|
| 352 | </p><p>
|
|---|
| 353 | </p></li><li><p>
|
|---|
| 354 | <span class="emphasis"><em>importance values > 0</em></span>
|
|---|
| 355 | </p><p>
|
|---|
| 356 | Normal allowed values
|
|---|
| 357 | </p><p>
|
|---|
| 358 | </p></li><li><p>
|
|---|
| 359 | <span class="emphasis"><em>Importance value smaller zero</em></span>
|
|---|
| 360 | </p><p>
|
|---|
| 361 | Not allowed!
|
|---|
| 362 | </p><p>
|
|---|
| 363 | </p></li></ul></div><p>
|
|---|
| 364 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.SampAlgor"></a>3.7.1.4.
|
|---|
| 365 | The Importance Sampling Algorithm
|
|---|
| 366 | </h4></div></div></div><p>
|
|---|
| 367 | Importance sampling supports using a customized importance
|
|---|
| 368 | sampling algorithm. To this end, the sampler interface
|
|---|
| 369 | <a href="ch03s07.html#anchor_EvtBias_G4VSampler">
|
|---|
| 370 | <code class="literal">G4VSampler</code></a>
|
|---|
| 371 | may be given a pointer to the interface
|
|---|
| 372 | <code class="literal">G4VImportanceAlgorithm</code>:
|
|---|
| 373 |
|
|---|
| 374 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 375 | class G4VImportanceAlgorithm
|
|---|
| 376 | {
|
|---|
| 377 | public:
|
|---|
| 378 | G4VImportanceAlgorithm();
|
|---|
| 379 | virtual ~G4VImportanceAlgorithm();
|
|---|
| 380 | virtual G4Nsplit_Weight Calculate(G4double ipre,
|
|---|
| 381 | G4double ipost,
|
|---|
| 382 | G4double init_w) const = 0;
|
|---|
| 383 | };
|
|---|
| 384 | </pre></div><p>
|
|---|
| 385 | </p><p>
|
|---|
| 386 | The method <code class="literal">Calculate()</code> takes the arguments:
|
|---|
| 387 |
|
|---|
| 388 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 389 | <span class="emphasis"><em>ipre, ipost</em></span>: importance
|
|---|
| 390 | of the previous cell and the importance of the current cell,
|
|---|
| 391 | respectively.
|
|---|
| 392 | </p></li><li><p>
|
|---|
| 393 | <span class="emphasis"><em>init_w</em></span>: the particles weight
|
|---|
| 394 | </p></li></ul></div><p>
|
|---|
| 395 | </p><p>
|
|---|
| 396 | It returns the struct:
|
|---|
| 397 |
|
|---|
| 398 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 399 | class G4Nsplit_Weight
|
|---|
| 400 | {
|
|---|
| 401 | public:
|
|---|
| 402 |
|
|---|
| 403 | G4int fN;
|
|---|
| 404 | G4double fW;
|
|---|
| 405 | };
|
|---|
| 406 | </pre></div><p>
|
|---|
| 407 |
|
|---|
| 408 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 409 | <span class="emphasis"><em>fN</em></span>: the calculated
|
|---|
| 410 | number of particles to exit the importance sampling
|
|---|
| 411 | </p></li><li><p>
|
|---|
| 412 | <span class="emphasis"><em>fW</em></span>: the weight of the particles
|
|---|
| 413 | </p></li></ul></div><p>
|
|---|
| 414 | </p><p>
|
|---|
| 415 | The user may have a customized algorithm used by providing a
|
|---|
| 416 | class inheriting from <code class="literal">G4VImportanceAlgorithm</code>.
|
|---|
| 417 | </p><p>
|
|---|
| 418 | If no customized algorithm is given to the sampler the default
|
|---|
| 419 | importance sampling algorithm is used. This algorithm is
|
|---|
| 420 | implemented in <code class="literal">G4ImportanceAlgorithm</code>.
|
|---|
| 421 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.WeightWin"></a>3.7.1.5.
|
|---|
| 422 | The Weight Window Technique
|
|---|
| 423 | </h4></div></div></div><p>
|
|---|
| 424 | The weight window technique is a weight-based alternative to
|
|---|
| 425 | importance sampling:
|
|---|
| 426 |
|
|---|
| 427 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 428 | applies splitting and Russian roulette depending on space
|
|---|
| 429 | (cells) and energy
|
|---|
| 430 | </p></li><li><p>
|
|---|
| 431 | user defines weight windows in contrast to defining importance
|
|---|
| 432 | values as in importance sampling
|
|---|
| 433 | </p></li></ul></div><p>
|
|---|
| 434 | </p><p>
|
|---|
| 435 | In contrast to importance sampling this technique is not weight
|
|---|
| 436 | blind. Instead the technique is applied according to the particle
|
|---|
| 437 | weight with respect to the current energy-space cell.
|
|---|
| 438 | </p><p>
|
|---|
| 439 | Therefore the technique is convenient to apply in combination
|
|---|
| 440 | with other variance reduction techniques such as cross-section
|
|---|
| 441 | biasing and implicit capture.
|
|---|
| 442 | </p><p>
|
|---|
| 443 | A weight window may be specified for every cell and for several
|
|---|
| 444 | energy regions: <span class="emphasis"><em>space-energy cell</em></span>.
|
|---|
| 445 |
|
|---|
| 446 | </p><div class="figure"><a name="fig.EvtBias.WeightWindow"></a><div class="figure-contents"><div class="mediaobject" align="center"><img src="./AllResources/Fundamentals/wwconcept.gif" align="middle" alt="Weight window concept"></div></div><p class="title"><b>Figure 3.2.
|
|---|
| 447 | Weight window concept
|
|---|
| 448 | </b></p></div><p><br class="figure-break">
|
|---|
| 449 | </p><h5><a name="id374547"></a>
|
|---|
| 450 | Weight window concept
|
|---|
| 451 | </h5><p>
|
|---|
| 452 | The user specifies a <span class="emphasis"><em>lower weight bound W_L</em></span>
|
|---|
| 453 | for every space-energy cell.
|
|---|
| 454 |
|
|---|
| 455 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 456 | The upper weight bound W_U and the survival weight W_S are
|
|---|
| 457 | calculated as:
|
|---|
| 458 | </p><p>
|
|---|
| 459 | W_U = C_U <span class="emphasis"><em>W_L</em></span> and
|
|---|
| 460 | </p><p>
|
|---|
| 461 | </p><p>
|
|---|
| 462 | W_S = C_S <span class="emphasis"><em>W_L</em></span>.
|
|---|
| 463 | </p><p>
|
|---|
| 464 | </p></li><li><p>
|
|---|
| 465 | The user specifies C_S and C_U once for the whole problem.
|
|---|
| 466 | </p></li><li><p>
|
|---|
| 467 | The user may give different sets of energy bounds for every cell
|
|---|
| 468 | or one set for all geometrical cells
|
|---|
| 469 | </p></li><li><p>
|
|---|
| 470 | Special case: if C_S = C_U = 1 for all energies then weight
|
|---|
| 471 | window is equivalent to importance sampling
|
|---|
| 472 | </p></li><li><p>
|
|---|
| 473 | The user can choose to apply the technique: at boundaries, on
|
|---|
| 474 | collisions or on boundaries and collisions
|
|---|
| 475 | </p></li></ul></div><p>
|
|---|
| 476 | </p><p>
|
|---|
| 477 | The energy-space cells are realized by <code class="literal">G4GeometryCell</code>
|
|---|
| 478 | as in importance sampling. The cells are stored in a weight window
|
|---|
| 479 | store defined by <code class="literal">G4VWeightWindowStore</code>:
|
|---|
| 480 |
|
|---|
| 481 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 482 | class G4VWeightWindowStore {
|
|---|
| 483 | public:
|
|---|
| 484 | G4VWeightWindowStore();
|
|---|
| 485 | virtual ~G4VWeightWindowStore();
|
|---|
| 486 | virtual G4double GetLowerWeitgh(const G4GeometryCell &gCell,
|
|---|
| 487 | G4double partEnergy) const = 0;
|
|---|
| 488 | virtual G4bool IsKnown(const G4GeometryCell &gCell) const = 0;
|
|---|
| 489 | virtual const G4VPhysicalVolume &GetWorldVolume() const = 0;
|
|---|
| 490 | };
|
|---|
| 491 | </pre></div><p>
|
|---|
| 492 | </p><p>
|
|---|
| 493 | A concrete implementation is provided:
|
|---|
| 494 |
|
|---|
| 495 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 496 | class G4WeightWindowStore: public G4VWeightWindowStore {
|
|---|
| 497 | public:
|
|---|
| 498 | explicit G4WeightWindowStore(const G4VPhysicalVolume &worldvolume);
|
|---|
| 499 | virtual ~G4WeightWindowStore();
|
|---|
| 500 | virtual G4double GetLowerWeitgh(const G4GeometryCell &gCell,
|
|---|
| 501 | G4double partEnergy) const;
|
|---|
| 502 | virtual G4bool IsKnown(const G4GeometryCell &gCell) const;
|
|---|
| 503 | virtual const G4VPhysicalVolume &GetWorldVolume() const;
|
|---|
| 504 | void AddLowerWeights(const G4GeometryCell &gCell,
|
|---|
| 505 | const std::vector<G4double> &lowerWeights);
|
|---|
| 506 | void AddUpperEboundLowerWeightPairs(const G4GeometryCell &gCell,
|
|---|
| 507 | const G4UpperEnergyToLowerWeightMap&
|
|---|
| 508 | enWeMap);
|
|---|
| 509 | void SetGeneralUpperEnergyBounds(const
|
|---|
| 510 | std::set<G4double, std::less<G4double> > & enBounds);
|
|---|
| 511 |
|
|---|
| 512 | private::
|
|---|
| 513 | ...
|
|---|
| 514 | };
|
|---|
| 515 | </pre></div><p>
|
|---|
| 516 | </p><p>
|
|---|
| 517 | The user may choose equal energy bounds for all cells. In this
|
|---|
| 518 | case a set of upper energy bounds must be given to the store using
|
|---|
| 519 | the method <code class="literal">SetGeneralUpperEnergyBounds</code>. If a general set
|
|---|
| 520 | of energy bounds have been set <code class="literal">AddLowerWeights</code> can be used
|
|---|
| 521 | to add the cells.
|
|---|
| 522 | </p><p>
|
|---|
| 523 | Alternatively, the user may chose different energy regions for
|
|---|
| 524 | different cells. In this case the user must provide a mapping of
|
|---|
| 525 | upper energy bounds to lower weight bounds for every cell using the
|
|---|
| 526 | method <code class="literal">AddUpperEboundLowerWeightPairs</code>.
|
|---|
| 527 | </p><p>
|
|---|
| 528 | Weight window algorithms implementing the interface class
|
|---|
| 529 | <code class="literal">G4VWeightWindowAlgorithm</code> can be used to define a
|
|---|
| 530 | customized algorithm:
|
|---|
| 531 |
|
|---|
| 532 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 533 | class G4VWeightWindowAlgorithm {
|
|---|
| 534 | public:
|
|---|
| 535 | G4VWeightWindowAlgorithm();
|
|---|
| 536 | virtual ~G4VWeightWindowAlgorithm();
|
|---|
| 537 | virtual G4Nsplit_Weight Calculate(G4double init_w,
|
|---|
| 538 | G4double lowerWeightBound) const = 0;
|
|---|
| 539 | };
|
|---|
| 540 | </pre></div><p>
|
|---|
| 541 | </p><p>
|
|---|
| 542 | A concrete implementation is provided and used as a default:
|
|---|
| 543 |
|
|---|
| 544 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 545 | class G4WeightWindowAlgorithm : public G4VWeightWindowAlgorithm {
|
|---|
| 546 | public:
|
|---|
| 547 | G4WeightWindowAlgorithm(G4double upperLimitFaktor = 5,
|
|---|
| 548 | G4double survivalFaktor = 3,
|
|---|
| 549 | G4int maxNumberOfSplits = 5);
|
|---|
| 550 | virtual ~G4WeightWindowAlgorithm();
|
|---|
| 551 | virtual G4Nsplit_Weight Calculate(G4double init_w,
|
|---|
| 552 | G4double lowerWeightBound) const;
|
|---|
| 553 | private:
|
|---|
| 554 | ...
|
|---|
| 555 | };
|
|---|
| 556 | </pre></div><p>
|
|---|
| 557 | </p><p>
|
|---|
| 558 | The constructor takes three parameters which are used to:
|
|---|
| 559 | calculate the upper weight bound (upperLimitFaktor), calculate the
|
|---|
| 560 | survival weight (survivalFaktor), and introduce a maximal number
|
|---|
| 561 | (maxNumberOfSplits) of copies to be created in one go.
|
|---|
| 562 | </p><p>
|
|---|
| 563 | In addition, the inverse of the maxNumberOfSplits is used to
|
|---|
| 564 | specify the minimum survival probability in case of Russian
|
|---|
| 565 | roulette.
|
|---|
| 566 | </p></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.ScorImpRoul.WeigRoul"></a>3.7.1.6.
|
|---|
| 567 | The Weight Roulette Technique
|
|---|
| 568 | </h4></div></div></div><p>
|
|---|
| 569 | Weight roulette (also called weight cutoff) is usually applied
|
|---|
| 570 | if importance sampling and implicit capture are used together.
|
|---|
| 571 | Implicit capture is not described here but it is useful to note
|
|---|
| 572 | that this procedure reduces a particle weight in every collision
|
|---|
| 573 | instead of killing the particle with some probability.
|
|---|
| 574 | </p><p>
|
|---|
| 575 | Together with importance sampling the weight of a particle may
|
|---|
| 576 | become so low that it does not change any result significantly.
|
|---|
| 577 | Hence tracking a very low weight particle is a waste of computing
|
|---|
| 578 | time. Weight roulette is applied in order to solve this
|
|---|
| 579 | problem.
|
|---|
| 580 | </p><h5><a name="id374768"></a>
|
|---|
| 581 | The weight roulette concept
|
|---|
| 582 | </h5><p>
|
|---|
| 583 | Weight roulette takes into account the importance "Ic" of the
|
|---|
| 584 | current cell and the importance "Is" of the cell in which the
|
|---|
| 585 | source is located, by using the ratio "R=Is/Ic".
|
|---|
| 586 | </p><p>
|
|---|
| 587 | Weight roulette uses a relative minimal weight limit and a
|
|---|
| 588 | relative survival weight. When a particle falls below the weight
|
|---|
| 589 | limit Russian roulette is applied. If the particle survives,
|
|---|
| 590 | tracking will be continued and the particle weight will be set to
|
|---|
| 591 | the survival weight.
|
|---|
| 592 | </p><p>
|
|---|
| 593 | The weight roulette uses the following parameters with their
|
|---|
| 594 | default values:
|
|---|
| 595 |
|
|---|
| 596 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 597 | <span class="emphasis"><em>wsurvival</em></span>: 0.5
|
|---|
| 598 | </p></li><li><p>
|
|---|
| 599 | <span class="emphasis"><em>wlimit</em></span>: 0.25
|
|---|
| 600 | </p></li><li><p>
|
|---|
| 601 | <span class="emphasis"><em>isource</em></span>: 1
|
|---|
| 602 | </p></li></ul></div><p>
|
|---|
| 603 | </p><p>
|
|---|
| 604 | The following algorithm is applied:
|
|---|
| 605 | </p><p>
|
|---|
| 606 | If a particle weight "w" is lower than R*wlimit:
|
|---|
| 607 |
|
|---|
| 608 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 609 | the weight of the particle will be changed to "ws = wsurvival*R"
|
|---|
| 610 | </p></li><li><p>
|
|---|
| 611 | the probability for the particle to survive is "p = w/ws"
|
|---|
| 612 | </p></li></ul></div><p>
|
|---|
| 613 | </p></div></div><div class="sect2" lang="en"><div class="titlepage"><div><div><h3 class="title"><a name="sect.EvtBias.PhysBias"></a>3.7.2.
|
|---|
| 614 | Physics Based Biasing
|
|---|
| 615 | </h3></div></div></div><p>
|
|---|
| 616 | Geant4 supports physics based biasing through a number of general
|
|---|
| 617 | use, built in biasing techniques. A utility class,
|
|---|
| 618 | G4WrapperProcess, is also available to support user defined
|
|---|
| 619 | biasing.
|
|---|
| 620 | </p><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.PhysBias.BuiltInOpt"></a>3.7.2.1.
|
|---|
| 621 | Built in Biasing Options
|
|---|
| 622 | </h4></div></div></div><div class="sect4" lang="en"><div class="titlepage"><div><div><h5 class="title"><a name="sect.EvtBias.PhysBias.BuiltInOpt.PrimPart"></a>3.7.2.1.1.
|
|---|
| 623 | Primary Particle Biasing
|
|---|
| 624 | </h5></div></div></div><p>
|
|---|
| 625 | Primary particle biasing can be used to increase the number of
|
|---|
| 626 | primary particles generated in a particular phase space region of
|
|---|
| 627 | interest. The weight of the primary particle is modified as
|
|---|
| 628 | appropriate. A general implementation is provided through the
|
|---|
| 629 | G4GeneralParticleSource class. It is possible to bias position,
|
|---|
| 630 | angular and energy distributions.
|
|---|
| 631 | </p><p>
|
|---|
| 632 | G4GeneralParticleSource is a concrete implementation of
|
|---|
| 633 | G4VPrimaryGenerator. To use, instantiate G4GeneralParticleSource in
|
|---|
| 634 | the G4VUserPrimaryGeneratorAction class, as demonstrated below.
|
|---|
| 635 |
|
|---|
| 636 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 637 | MyPrimaryGeneratorAction::MyPrimaryGeneratorAction() {
|
|---|
| 638 | generator = new G4GeneralParticleSource;
|
|---|
| 639 | }
|
|---|
| 640 |
|
|---|
| 641 | void
|
|---|
| 642 | MyPrimaryGeneratorAction::GeneratePrimaries(G4Event*anEvent){
|
|---|
| 643 | generator->GeneratePrimaryVertex(anEvent);
|
|---|
| 644 | }
|
|---|
| 645 | </pre></div><p>
|
|---|
| 646 | </p><p>
|
|---|
| 647 | The biasing can be configured through interactive commands.
|
|---|
| 648 | Extensive documentation can be found in
|
|---|
| 649 | <a href="http://reat.space.qinetiq.com/gps/" target="_top">
|
|---|
| 650 | Primary particle biasing</a>. Examples are also distributed
|
|---|
| 651 | with the Geant4 distribution in
|
|---|
| 652 | <span class="bold"><strong>examples/extended/eventgenerator/exgps</strong></span>.
|
|---|
| 653 | </p></div><div class="sect4" lang="en"><div class="titlepage"><div><div><h5 class="title"><a name="sect.EvtBias.PhysBias.BuiltInOpt.RadDcy"></a>3.7.2.1.2.
|
|---|
| 654 | Radioactive Decay Biasing
|
|---|
| 655 | </h5></div></div></div><p>
|
|---|
| 656 | The G4RadioactiveDecay class simulates the decay of radioactive
|
|---|
| 657 | nuclei and implements the following biasing options:
|
|---|
| 658 |
|
|---|
| 659 | </p><div class="itemizedlist"><ul type="disc" compact><li><p>
|
|---|
| 660 | Increase the sampling rate of radionuclides within observation
|
|---|
| 661 | times through a user defined probability distribution function
|
|---|
| 662 | </p></li><li><p>
|
|---|
| 663 | Nuclear splitting, where the parent nuclide is split into a
|
|---|
| 664 | user defined number of nuclides
|
|---|
| 665 | </p></li><li><p>
|
|---|
| 666 | Branching ratio biasing where branching ratios are sampled with
|
|---|
| 667 | equal probability
|
|---|
| 668 | </p></li></ul></div><p>
|
|---|
| 669 | </p><p>
|
|---|
| 670 | G4RadioactiveDecay is a process which must be registered with a
|
|---|
| 671 | process manager, as demonstrated below.
|
|---|
| 672 |
|
|---|
| 673 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 674 | void MyPhysicsList::ConstructProcess()
|
|---|
| 675 | {
|
|---|
| 676 | ...
|
|---|
| 677 | G4RadioactiveDecay* theRadioactiveDecay =
|
|---|
| 678 | new G4RadioactiveDecay();
|
|---|
| 679 |
|
|---|
| 680 | G4ProcessManager* pmanager = ...
|
|---|
| 681 | pmanager ->AddProcess(theRadioactiveDecay);
|
|---|
| 682 | ...
|
|---|
| 683 | }
|
|---|
| 684 | </pre></div><p>
|
|---|
| 685 | </p><p>
|
|---|
| 686 | The biasing can be controlled either in compiled code or through
|
|---|
| 687 | interactive commands. Extensive documentation can be found in
|
|---|
| 688 |
|
|---|
| 689 | <a href="http://reat.space.qinetiq.com/septimess/exrdm/" target="_top">
|
|---|
| 690 | Radioactive decay biasing example
|
|---|
| 691 | </a>
|
|---|
| 692 | and
|
|---|
| 693 | <a href="http://www.space.qinetiq.com/geant4/rdm.html" target="_top">
|
|---|
| 694 | Radioactive decay biasing
|
|---|
| 695 | </a>.
|
|---|
| 696 | </p><p>
|
|---|
| 697 | Radioactive decay biasing examples are also distributed with the Geant4
|
|---|
| 698 | distribution in
|
|---|
| 699 | <span class="bold"><strong>examples/extended/radioactivedecay/exrdm</strong></span>.
|
|---|
| 700 | </p></div><div class="sect4" lang="en"><div class="titlepage"><div><div><h5 class="title"><a name="sect.EvtBias.PhysBias.BuiltInOpt.HadLeadPar"></a>3.7.2.1.3.
|
|---|
| 701 | Hadronic Leading Particle Biasing
|
|---|
| 702 | </h5></div></div></div><p>
|
|---|
| 703 | One hadronic leading particle biasing technique is
|
|---|
| 704 | implemented in the G4HadLeadBias utility. This method keeps only
|
|---|
| 705 | the most important part of the event, as well as representative
|
|---|
| 706 | tracks of each given particle type. So the track with the highest
|
|---|
| 707 | energy as well as one of each of Baryon, pi0, mesons and leptons.
|
|---|
| 708 | As usual, appropriate weights are assigned to the particles.
|
|---|
| 709 | Setting the <span class="bold"><strong>SwitchLeadBiasOn</strong></span>
|
|---|
| 710 | environmental variable will activate this utility.
|
|---|
| 711 | </p></div><div class="sect4" lang="en"><div class="titlepage"><div><div><h5 class="title"><a name="sect.EvtBias.PhysBias.BuiltInOpt.HadCrsSect"></a>3.7.2.1.4.
|
|---|
| 712 | Hadronic Cross Section Biasing
|
|---|
| 713 | </h5></div></div></div><p>
|
|---|
| 714 | Cross section biasing artificially enhances/reduces the cross
|
|---|
| 715 | section of a process. This may be useful for studying thin layer
|
|---|
| 716 | interactions or thick layer shielding. The built in hadronic cross
|
|---|
| 717 | section biasing applies to photon inelastic, electron nuclear and
|
|---|
| 718 | positron nuclear processes.
|
|---|
| 719 | </p><p>
|
|---|
| 720 | The biasing is controlled through the
|
|---|
| 721 | <span class="bold"><strong>BiasCrossSectionByFactor</strong></span> method
|
|---|
| 722 | in G4HadronicProcess, as demonstrated below.
|
|---|
| 723 |
|
|---|
| 724 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 725 | void MyPhysicsList::ConstructProcess()
|
|---|
| 726 | {
|
|---|
| 727 | ...
|
|---|
| 728 | G4ElectroNuclearReaction * theElectroReaction =
|
|---|
| 729 | new G4ElectroNuclearReaction;
|
|---|
| 730 |
|
|---|
| 731 | G4ElectronNuclearProcess theElectronNuclearProcess;
|
|---|
| 732 | theElectronNuclearProcess.RegisterMe(theElectroReaction);
|
|---|
| 733 | theElectronNuclearProcess.BiasCrossSectionByFactor(100);
|
|---|
| 734 |
|
|---|
| 735 | pManager->AddDiscreteProcess(&theElectronNuclearProcess);
|
|---|
| 736 | ...
|
|---|
| 737 | }
|
|---|
| 738 | </pre></div><p>
|
|---|
| 739 | </p><p>
|
|---|
| 740 | More details can be found in
|
|---|
| 741 | <a href="http://www.triumf.ca/geant4-03/talks/03-Wednesday-AM-1/03-J.Wellisch/biasing.hadronics.pdf" target="_top">
|
|---|
| 742 | Hadronic cross section biasing
|
|---|
| 743 | </a>.
|
|---|
| 744 | </p></div></div><div class="sect3" lang="en"><div class="titlepage"><div><div><h4 class="title"><a name="sect.EvtBias.PhysBias."></a>3.7.2.2.
|
|---|
| 745 | G4WrapperProcess
|
|---|
| 746 | </h4></div></div></div><p>
|
|---|
| 747 | G4WrapperProcess can be used to implement user defined event
|
|---|
| 748 | biasing. G4WrapperProcess, which is a process itself, wraps an
|
|---|
| 749 | existing process. By default, all function calls are forwared to
|
|---|
| 750 | the wrapped process. It is a non-invasive way to modify the
|
|---|
| 751 | behaviour of an existing process.
|
|---|
| 752 | </p><p>
|
|---|
| 753 | To use this utility, first create a derived class inheriting
|
|---|
| 754 | from G4WrapperProcess. Override the methods whose behaviour you
|
|---|
| 755 | would like to modify, for example, PostStepDoIt, and register the
|
|---|
| 756 | derived class in place of the process to be wrapped. Finally,
|
|---|
| 757 | register the wrapped process with G4WrapperProcess. The code
|
|---|
| 758 | snippets below demonstrate its use.
|
|---|
| 759 |
|
|---|
| 760 | </p><div class="informalexample"><pre class="programlisting">
|
|---|
| 761 | class MyWrapperProcess : public G4WrapperProcess {
|
|---|
| 762 | ...
|
|---|
| 763 | G4VParticleChange* PostStepDoIt(const G4Track& track,
|
|---|
| 764 | const G4Step& step) {
|
|---|
| 765 | // Do something interesting
|
|---|
| 766 | }
|
|---|
| 767 | };
|
|---|
| 768 |
|
|---|
| 769 |
|
|---|
| 770 | void MyPhysicsList::ConstructProcess()
|
|---|
| 771 | {
|
|---|
| 772 | ...
|
|---|
| 773 | G4LowEnergyBremsstrahlung* bremProcess =
|
|---|
| 774 | new G4LowEnergyBremsstrahlung();
|
|---|
| 775 |
|
|---|
| 776 | MyWrapperProcess* wrapper = new MyWrapperProcess();
|
|---|
| 777 | wrapper->RegisterProcess(bremProcess);
|
|---|
| 778 |
|
|---|
| 779 | processManager->AddProcess(wrapper, -1, -1, 3);
|
|---|
| 780 | }
|
|---|
| 781 | </pre></div><p>
|
|---|
| 782 | </p></div></div></div><div class="navfooter"><hr><table width="100%" summary="Navigation footer"><tr><td width="40%" align="left"><a accesskey="p" href="ch03s06.html"><img src="AllResources/IconsGIF/prev.gif" alt="Prev"></a> </td><td width="20%" align="center"><a accesskey="u" href="ch03.html"><img src="AllResources/IconsGIF/up.gif" alt="Up"></a></td><td width="40%" align="right"> <a accesskey="n" href="ch04.html"><img src="AllResources/IconsGIF/next.gif" alt="Next"></a></td></tr><tr><td width="40%" align="left" valign="top">3.6.
|
|---|
| 783 | Event Generator Interface
|
|---|
| 784 | </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"> Chapter 4.
|
|---|
| 785 | Detector Definition and Response
|
|---|
| 786 | </td></tr></table></div></body></html>
|
|---|