Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/biasing/GB01/include/GB01BOptrMultiParticleChangeCrossSection.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 GB01/include/GB01BOptrMultiParticleChangeCrossSection.hh
 27 /// \brief Definition of the GB01BOptrMultiParticleChangeCrossSection class
 28 //
 29 //-----------------------------------------------------------------
 30 //
 31 // GB01BOptrMultiParticleChangeCrossSection
 32 //
 33 // Class Description:
 34 //        A G4VBiasingOperator concrete implementation example to
 35 //    illustrate how an operator can make use of an other operator.
 36 //        In the present case, the GB01BOptrChangeCrossSection
 37 //    operator, that is made to bias processes of one particle
 38 //    type is instancied several times, one for each particle type
 39 //    to bias.
 40 //
 41 //-----------------------------------------------------------------
 42 //
 43 
 44 #ifndef GB01BOptrMultiParticleChangeCrossSection_hh
 45 #define GB01BOptrMultiParticleChangeCrossSection_hh 1
 46 
 47 #include "G4VBiasingOperator.hh"
 48 class GB01BOptrChangeCrossSection;
 49 class G4ParticleDefinition;
 50 
 51 #include <map>
 52 
 53 class GB01BOptrMultiParticleChangeCrossSection : public G4VBiasingOperator
 54 {
 55   public:
 56     GB01BOptrMultiParticleChangeCrossSection();
 57     virtual ~GB01BOptrMultiParticleChangeCrossSection() {}
 58 
 59     // ---------------------------------
 60     // -- Method specific to this class:
 61     // ---------------------------------
 62     // -- Each particle type for which its name is passed will be biased; *provided*
 63     // -- that the proper calls to biasingPhysics->Bias(particleName) have been done
 64     // -- in the main program.
 65     void AddParticle(G4String particleName);
 66 
 67   private:
 68     // -----------------------------
 69     // -- Mandatory from base class:
 70     // -----------------------------
 71     // -- This method returns a biasing operation that will bias the physics process occurence:
 72     virtual G4VBiasingOperation*
 73     ProposeOccurenceBiasingOperation(const G4Track* track,
 74                                      const G4BiasingProcessInterface* callingProcess);
 75     // -- Methods not used:
 76     virtual G4VBiasingOperation* ProposeFinalStateBiasingOperation(const G4Track*,
 77                                                                    const G4BiasingProcessInterface*)
 78     {
 79       return 0;
 80     }
 81     virtual G4VBiasingOperation* ProposeNonPhysicsBiasingOperation(const G4Track*,
 82                                                                    const G4BiasingProcessInterface*)
 83     {
 84       return 0;
 85     }
 86 
 87   private:
 88     // -- ("using" is to avoid compiler complaining against (false) method shadowing.)
 89     using G4VBiasingOperator::OperationApplied;
 90 
 91     // -- Optionnal base class method implementation.
 92     // -- This method is called to inform the operator that a proposed operation has been applied.
 93     // -- In the present case, it means that a physical interaction occured (interaction at
 94     // -- PostStepDoIt level):
 95     virtual void OperationApplied(const G4BiasingProcessInterface* callingProcess,
 96                                   G4BiasingAppliedCase biasingCase,
 97                                   G4VBiasingOperation* occurenceOperationApplied,
 98                                   G4double weightForOccurenceInteraction,
 99                                   G4VBiasingOperation* finalStateOperationApplied,
100                                   const G4VParticleChange* particleChangeProduced);
101 
102   public:
103     // -- Optionnal base class method. It is called at the time a tracking of a particle starts:
104     void StartTracking(const G4Track* track);
105 
106   private:
107     // -- List of associations between particle types and biasing operators:
108     std::map<const G4ParticleDefinition*, GB01BOptrChangeCrossSection*> fBOptrForParticle;
109     std::vector<const G4ParticleDefinition*> fParticlesToBias;
110     GB01BOptrChangeCrossSection* fCurrentOperator;
111 
112     // -- count number of biased interations for current track:
113     G4int fnInteractions;
114 };
115 
116 #endif
117