Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/biasing/GB06/include/GB06BOptrSplitAndKillByImportance.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 /// \file GB06/include/GB06BOptrSplitAndKillByImportance.hh
 27 /// \brief Definition of the GB06BOptrSplitAndKillByImportance class
 28 //
 29 //---------------------------------------------------------------
 30 //
 31 // GB06BOptrSplitAndKillByImportance
 32 //
 33 // Class Description:
 34 //        A G4VBiasingOperator concrete implementation example to
 35 //    illustrate how to bias physics processes cross-section for
 36 //    one particle type.
 37 //        The G4VBiasingOperation G4BOptnChangeImportance is
 38 //    selected by this operator, and is sent to each process
 39 //    calling the operator.
 40 //        A simple constant bias to the cross-section is applied,
 41 //    but more sophisticated changes can be applied.
 42 //
 43 //---------------------------------------------------------------
 44 //
 45 
 46 #ifndef GB06BOptrSplitAndKillByImportance_hh
 47 #define GB06BOptrSplitAndKillByImportance_hh 1
 48 
 49 #include "G4VBiasingOperator.hh"
 50 
 51 class GB06BOptnSplitAndKillByImportance;
 52 class G4BiasingProcessSharedData;
 53 class G4ParallelGeometriesLimiterProcess;
 54 class G4ParticleDefinition;
 55 class G4VPhysicalVolume;
 56 #include <map>
 57 
 58 class GB06BOptrSplitAndKillByImportance : public G4VBiasingOperator
 59 {
 60   public:
 61     // ------------------------------------------------------------
 62     // -- Constructor: takes the name of the particle type to bias:
 63     // ------------------------------------------------------------
 64     GB06BOptrSplitAndKillByImportance(G4String particleToBias,
 65                                       G4String name = "SplitAndKillByImportance");
 66     virtual ~GB06BOptrSplitAndKillByImportance();
 67 
 68     // -- method called at beginning of run:
 69     virtual void StartRun();
 70 
 71   private:
 72     // -----------------------------
 73     // -- Mandatory from base class:
 74     // -----------------------------
 75     // -- Not used:
 76     virtual G4VBiasingOperation*
 77     ProposeOccurenceBiasingOperation(const G4Track*, const G4BiasingProcessInterface*) final
 78     {
 79       return 0;
 80     }
 81 
 82     // -- Not used:
 83     virtual G4VBiasingOperation*
 84     ProposeFinalStateBiasingOperation(const G4Track*, const G4BiasingProcessInterface*) final
 85     {
 86       return 0;
 87     }
 88 
 89     // -- Used method : it will return the biasing operation that will split particles
 90     // -- with a probabilty depending on the total absorption cross-section.
 91     virtual G4VBiasingOperation*
 92     ProposeNonPhysicsBiasingOperation(const G4Track*, const G4BiasingProcessInterface*) final;
 93 
 94     // ---------------------------------------
 95     // -- Method specific to this application:
 96     // ---------------------------------------
 97   public:
 98     void SetParallelWorld(G4VPhysicalVolume* parallelWorld) { fParallelWorld = parallelWorld; }
 99 
100     // -- The importance map, linking a replica number to an volume importance:
101     std::map<G4int, G4int>& GetImportanceMap() { return fImportanceMap; }
102 
103   private:
104     GB06BOptnSplitAndKillByImportance* fSplitAndKillByImportance;
105     const G4ParticleDefinition* fParticleToBias;
106     G4VPhysicalVolume* fParallelWorld;
107     G4int fParallelWorldIndex;
108     const G4BiasingProcessSharedData* fBiasingSharedData;
109     const G4ParallelGeometriesLimiterProcess* fBiasingLimiterProcess;
110     std::map<G4int, G4int> fImportanceMap;
111 };
112 
113 #endif
114