ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4IT.hh
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4IT.hh
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 // Author: Mathieu Karamitros
28 
29 // The code is developed in the framework of the ESA AO7146
30 //
31 // We would be very happy hearing from you, send us your feedback! :)
32 //
33 // In order for Geant4-DNA to be maintained and still open-source,
34 // article citations are crucial.
35 // If you use Geant4-DNA chemistry and you publish papers about your software,
36 // in addition to the general paper on Geant4-DNA:
37 //
38 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39 //
40 // we would be very happy if you could please also cite the following
41 // reference papers on chemistry:
42 //
43 // J. Comput. Phys. 274 (2014) 841-882
44 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45 
46 #ifndef G4IT_h
47 #define G4IT_h 1
48 
49 #include "globals.hh"
50 #include "G4ITType.hh"
51 #include "G4ThreeVector.hh"
53 #include "G4KDNode.hh"
54 
56 // To implement your own IT class, you should use
57 // ITDef(MyClass) in the class define in your MyClass.hh
58 // and ITImp(MyClass) in your MyClass.cc
59 // For instance, see G4Molecule
61 
62 class G4IT;
63 template<>
65 
67 //template<typename PointT> class G4KDNode;
68 class G4KDNode_Base;
69 class G4ITBox;
70 class G4Track;
71 
72 G4IT* GetIT(const G4Track* track);
73 G4IT* GetIT(const G4Track& track);
74 
75 template<class OBJECT> class G4FastListNode;
77 
87 class G4IT : public virtual G4VUserTrackInformation
88 {
89 public:
90  G4IT();
91  G4IT(G4Track*);
92  virtual ~G4IT();
93 
94 // inline void *operator new(size_t);
95 // inline void operator delete(void *aIT);
96 
97  virtual void Print() const
98  {
99  ;
100  }
101  virtual const G4String& GetName() const = 0;
102 
104  // You should not worried of implementing diff, equal
105  // and GetType.
106  // When using ITDef(MyClass) this will be done.
107  // However, you need to implement in the concrete class
108  // even fake operators for < and ==
109  // They will be used by diff and equal.
111  virtual G4bool diff(const G4IT& right) const = 0;
112  virtual G4bool equal(const G4IT& right) const = 0;
113  G4bool operator<(const G4IT& right) const;
114  G4bool operator==(const G4IT& right) const;
115  G4bool operator!=(const G4IT& right) const;
116 
117  void SetTrack(G4Track*);
118  inline G4Track* GetTrack();
119  inline const G4Track* GetTrack() const;
120 
122  const G4ThreeVector& GetPosition() const;
123  double operator[](int i) const;
124  const G4ThreeVector& GetPreStepPosition() const;
127 
128  inline void SetPrevious(G4IT*);
129  inline void SetNext(G4IT*);
130  inline G4IT* GetPrevious();
131  inline G4IT* GetNext();
132  inline const G4IT* GetPrevious() const;
133  inline const G4IT* GetNext() const;
134  inline void SetITBox(G4ITBox*);
135  inline const G4ITBox* GetITBox() const;
136  void TakeOutBox();
137  inline void SetNode(G4KDNode_Base*);
138  inline G4KDNode_Base* GetNode() const;
139 
140  inline void SetParentID(int, int);
141  inline void GetParentID(int&, int&);
142 
144  {
145  return fpTrackingInformation;
146  }
147 
149  {
150  return fpTrackNode;
151  }
152  inline void SetListNode(G4TrackListNode* node)
153  {
154  fpTrackNode = node;
155  }
156 
157  virtual const G4ITType GetITType() const = 0;
158 
159  virtual G4ITType GetITSubType() const
160  {
161  return 0;
162  }
163 
164 protected:
165  G4IT(const G4IT&);
166  G4IT& operator=(const G4IT&);
168 
169 private:
174 
177 
180 };
181 //------------------------------------------------------------------------------
182 
183 inline const G4ITBox* G4IT::GetITBox() const
184 {
185  return fpITBox;
186 }
187 
188 inline void G4IT::SetITBox(G4ITBox * aITBox)
189 {
190  fpITBox = aITBox;
191 }
192 
193 inline void G4IT::SetPrevious(G4IT* aIT)
194 {
195  fpPreviousIT = aIT;
196 }
197 
198 inline void G4IT::SetNext(G4IT* aIT)
199 {
200  fpNextIT = aIT;
201 }
202 
204 {
205  return fpPreviousIT;
206 }
207 
209 {
210  return fpNextIT;
211 }
212 
214 {
215  fpTrack = track;
216 }
217 
219 {
220  return fpTrack;
221 }
222 
223 inline const G4Track* G4IT::GetTrack() const
224 {
225  return fpTrack;
226 }
227 
228 inline void G4IT::SetParentID(int p_a, int p_b)
229 {
230  fParentID_A = p_a;
231  fParentID_B = p_b;
232 }
233 
234 inline void G4IT::GetParentID(int& p_a, int&p_b)
235 {
236  p_a = fParentID_A;
237  p_b = fParentID_B;
238 }
239 
240 inline const G4IT* G4IT::GetPrevious() const
241 {
242  return fpPreviousIT;
243 }
244 
245 inline const G4IT* G4IT::GetNext() const
246 {
247  return fpNextIT;
248 }
249 
250 inline void G4IT::SetNode(G4KDNode_Base* aNode)
251 {
252  fpKDNode = aNode;
253 }
254 
256 {
257  return fpKDNode;
258 }
259 #endif