Geant4  10.03
G4WrapperProcess.hh
Go to the documentation of this file.
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 //
27 // $Id: G4WrapperProcess.hh 80787 2014-05-12 09:06:07Z gcosmo $
28 //
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class header file
32 //
33 // Class Description
34 //
35 // This class is the virtual class for wrapper process objects.
36 
37 // ------------------------------------------------------------
38 // New Physics scheme 18 Dec. 1996 H.Kurahige
39 // ------------------------------------------------------------
40 
41 #ifndef G4WrapperProcess_h
42 #define G4WrapperProcess_h 1
43 
44 #include "globals.hh"
45 #include "G4ios.hh"
46 #include "G4VProcess.hh"
47 
49 {
50  // A virtual class for wrapper process objects.
51 
52  private:
53  // hide default constructor and assignment operator as private
54  // do not hide default constructor for alpha version
56 
57  public: // with description
58  // constructor requires the process name and type
59  G4WrapperProcess(const G4String& aName = "Wrapped",
60  G4ProcessType aType = fNotDefined );
61 
62  // copy constructor copys the name but does not copy the
63  // physics table (0 pointer is assigned)
64  G4WrapperProcess(const G4WrapperProcess &right);
65 
66  public:
67  // destructor
68  virtual ~G4WrapperProcess();
69 
70  // equality opperators
71  inline G4int operator==(const G4WrapperProcess &right) const;
72  inline G4int operator!=(const G4WrapperProcess &right) const;
73 
74  public: // with description
75  virtual void RegisterProcess(G4VProcess*);
76  virtual const G4VProcess* GetRegisteredProcess() const;
77 
78  protected:
80 
81  public: // with description
83  // DoIt /////////////////
86  const G4Track& track,
87  const G4Step& stepData
88  );
89 
91  const G4Track& track,
92  const G4Step& stepData
93  );
95  const G4Track& track,
96  const G4Step& stepData
97  );
99  // GPIL //////////////
102  const G4Track& track,
103  G4double previousStepSize,
104  G4double currentMinimumStep,
105  G4double& proposedSafety,
106  G4GPILSelection* selection);
107 
109  const G4Track& track,
111  );
112 
114  const G4Track& track,
115  G4double previousStepSize,
117  ) ;
118 
120  virtual G4bool IsApplicable(const G4ParticleDefinition&);
121  // Returns true if this process object is applicable to
122  // the particle type
123  // Process will not be registered to a particle if IsApplicable is false
124 
125  virtual void BuildPhysicsTable(const G4ParticleDefinition&);
126  // Messaged by the Particle definition (via the Process manager)
127  // whenever cross section tables have to be rebuilt (i.e. if new
128  // materials have been defined).
129  // It is overloaded by individual processes when they need physics
130  // tables.
131 
132  // Processes which Build (for example in their
133  // constructors) physics tables independent of cuts
134  // should preferably use a
135  // private void BuildThePhysicsTable()
136  // function. Not another BuildPhysicsTable, please.
137 
138  virtual void PreparePhysicsTable(const G4ParticleDefinition&);
139  // Messaged by the Particle definition (via the Process manager)
140  // whenever cross section tables have to be prepare for rebuilt
141  // (i.e. if new materials have been defined).
142  // It is overloaded by individual processes when they need physics
143  // tables.
144 
145  // Processes which Build physics tables independent of cuts
146  // (for example in their constructors)
147  // should preferably use private
148  // void BuildThePhysicsTable() and void PreparePhysicsTable().
149  // Not another BuildPhysicsTable, please.
150 
151 
153  const G4String& directory,
154  G4bool ascii = false);
155  // Store PhysicsTable in a file.
156  // (return false in case of failure at I/O )
157 
159  const G4String& directory,
160  G4bool ascii = false);
161  // Retrieve Physics from a file.
162  // (return true if the Physics Table can be build by using file)
163  // (return false if the process has no functionality or in case of failure)
164  // File name should be defined by each process
165  // and the file should be placed under the directory specifed by the argument.
167  virtual void StartTracking(G4Track*);
168  virtual void EndTracking();
169  // inform Start/End of tracking for each track to the physics process
170 
171  public:
172  virtual void SetProcessManager(const G4ProcessManager*);
173  // A process manager set its own pointer when the process is registered
174  // the process Manager
175  virtual const G4ProcessManager* GetProcessManager();
176  // Get the process manager which the process belongs to
177 
178  public:
179  virtual void ResetNumberOfInteractionLengthLeft();
180  // reset (determine the value of)NumberOfInteractionLengthLeft
181  virtual void SetMasterProcess(G4VProcess* masterP);
182  // Needed for MT, forward call to underlying process
183 };
184 
185 inline
187 {
188  G4Exception("G4WrapperProcess::operator=","Illegal operation",
189  JustWarning,"Assignment operator is called");
190  return *this;
191 }
192 
193 inline
195 {
196  return (this == &right);
197 }
198 
199 inline
201 {
202  return (this != &right);
203 }
204 
205 #endif
G4int operator!=(const G4WrapperProcess &right) const
G4double condition(const G4ErrorSymMatrix &m)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)
virtual void RegisterProcess(G4VProcess *)
G4WrapperProcess(const G4String &aName="Wrapped", G4ProcessType aType=fNotDefined)
virtual void ResetNumberOfInteractionLengthLeft()
int G4int
Definition: G4Types.hh:78
virtual void SetProcessManager(const G4ProcessManager *)
virtual G4bool IsApplicable(const G4ParticleDefinition &)
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
virtual ~G4WrapperProcess()
virtual G4VParticleChange * AtRestDoIt(const G4Track &track, const G4Step &stepData)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
G4int operator==(const G4WrapperProcess &right) const
G4VProcess * pRegProcess
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual const G4VProcess * GetRegisteredProcess() const
virtual void EndTracking()
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
virtual void SetMasterProcess(G4VProcess *masterP)
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
double G4double
Definition: G4Types.hh:76
G4ForceCondition
virtual const G4ProcessManager * GetProcessManager()
virtual void StartTracking(G4Track *)
G4WrapperProcess & operator=(const G4WrapperProcess &right)
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4GPILSelection
G4ProcessType