Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PrimaryVertex.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// G4PrimaryVertex
27//
28// class description:
29//
30// This is the class which represents a primary vertex. The object of this
31// class is set to G4Event by a G4VPrimaryGenerator concrete class.
32// This class object has one or more G4PrimaryParticle objects as primary
33// particles.
34
35// Authors: G.Cosmo, 2 December 1995 - Design, based on object model
36// M.Asai, 29 January 1996 - First implementation
37// --------------------------------------------------------------------
38#ifndef G4PrimaryVertex_hh
39#define G4PrimaryVertex_hh 1
40
41#include "globals.hh"
42#include "pwdefs.hh"
43#include "G4Allocator.hh"
44#include "G4ThreeVector.hh"
45#include "G4PrimaryParticle.hh"
46
48
50{
51 public:
52
56 // Constructors
57
58 virtual ~G4PrimaryVertex();
59 // Destructor
60
61 G4PrimaryVertex(const G4PrimaryVertex& right);
63 // Copy constructor & assignment operator
64
65 G4bool operator==(const G4PrimaryVertex& right) const;
66 G4bool operator!=(const G4PrimaryVertex& right) const;
67 // Equality operators
68
69 inline void* operator new(size_t);
70 inline void operator delete(void* aPrimaryVertex);
71 // Overloaded new/delete operators
72
73 // Accessors & modifiers
74
75 inline G4ThreeVector GetPosition() const;
76 inline void SetPosition(G4double x0, G4double y0, G4double z0);
77 inline G4double GetX0() const;
78 inline G4double GetY0() const;
79 inline G4double GetZ0() const;
80 inline G4double GetT0() const;
81 inline void SetT0(G4double t0);
82 inline G4int GetNumberOfParticle() const;
83 inline void SetPrimary(G4PrimaryParticle* pp);
85 inline void SetNext(G4PrimaryVertex* nv);
86 inline void ClearNext();
87 inline G4PrimaryVertex* GetNext() const;
88 inline G4double GetWeight() const;
89 inline void SetWeight(G4double w);
92
93 void Print() const;
94
95 private:
96
97 G4double X0 = 0.0;
98 G4double Y0 = 0.0;
99 G4double Z0 = 0.0;
100 G4double T0 = 0.0;
101 G4PrimaryParticle* theParticle = nullptr;
102 G4PrimaryParticle* theTail = nullptr;
103 G4PrimaryVertex* nextVertex = nullptr;
104 G4PrimaryVertex* tailVertex = nullptr;
105 G4double Weight0 = 1.0;
106 G4VUserPrimaryVertexInformation* userInfo = nullptr;
107 G4int numberOfParticle = 0;
108};
109
111
112// ------------------------
113// Inline methods
114// ------------------------
115
116inline
117void* G4PrimaryVertex::operator new(std::size_t)
118{
119 if (aPrimaryVertexAllocator() == nullptr)
120 {
122 }
123 return (void*) aPrimaryVertexAllocator()->MallocSingle();
124}
125
126inline
127void G4PrimaryVertex::operator delete(void* aPrimaryVertex)
128{
129 aPrimaryVertexAllocator()->FreeSingle((G4PrimaryVertex*) aPrimaryVertex);
130}
131
132inline
134{
135 return G4ThreeVector(X0,Y0,Z0);
136}
137
138inline
140{
141 X0 = x0; Y0 = y0; Z0 = z0;
142}
143
144inline
146{
147 return X0;
148}
149
150inline
152{
153 return Y0;
154}
155
156inline
158{
159 return Z0;
160}
161
162inline
164{
165 return T0;
166}
167
168inline
170{
171 T0 = t0;
172}
173
174inline
176{
177 return numberOfParticle;
178}
179
180inline
182{
183 if(theParticle == nullptr) { theParticle = pp; }
184 else { theTail->SetNext(pp); }
185 theTail = pp;
186 ++numberOfParticle;
187}
188
189inline
191{
192 if(nextVertex == nullptr) { nextVertex = nv; }
193 else { tailVertex->SetNext(nv); }
194 tailVertex = nv;
195}
196
197inline
199{
200 nextVertex = nullptr;
201 tailVertex = nullptr;
202}
203
204inline
206{
207 return nextVertex;
208}
209
210inline
212{
213 return Weight0;
214}
215
216inline
218{
219 Weight0 = w;
220}
221
222inline
224{
225 userInfo = info;
226}
227
228inline
230{
231 return userInfo;
232}
233
234#endif
235
G4PART_DLL G4Allocator< G4PrimaryVertex > *& aPrimaryVertexAllocator()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetNext(G4PrimaryParticle *np)
G4double GetT0() const
G4bool operator==(const G4PrimaryVertex &right) const
void Print() const
virtual ~G4PrimaryVertex()
void SetNext(G4PrimaryVertex *nv)
G4PrimaryVertex * GetNext() const
G4double GetWeight() const
G4VUserPrimaryVertexInformation * GetUserInformation() const
G4PrimaryVertex & operator=(const G4PrimaryVertex &right)
G4double GetZ0() const
G4bool operator!=(const G4PrimaryVertex &right) const
void SetPosition(G4double x0, G4double y0, G4double z0)
G4double GetX0() const
G4double GetY0() const
G4ThreeVector GetPosition() const
void SetPrimary(G4PrimaryParticle *pp)
void SetWeight(G4double w)
void SetT0(G4double t0)
void SetUserInformation(G4VUserPrimaryVertexInformation *info)
G4PrimaryParticle * GetPrimary(G4int i=0) const
G4int GetNumberOfParticle() const
#define G4PART_DLL
Definition: pwdefs.hh:45