source: trunk/source/processes/hadronic/models/util/include/G4Fragment.hh

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

File size: 10.0 KB
RevLine 
[819]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//
[1340]26// $Id: G4Fragment.hh,v 1.16 2010/09/28 16:09:00 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-ref-09 $
[819]28//
[1315]29//---------------------------------------------------------------------
[819]30//
[1315]31// Geant4 header G4Fragment
32//
[819]33// by V. Lara (May 1998)
[1315]34//
35// Modifications:
36// 03.05.2010 V.Ivanchenko General cleanup of inline functions: objects
37//            are accessed by reference; remove double return
38//            tolerance of excitation energy at modent it is computed;
39//            safe computation of excitation for exotic fragments
40// 18.05.2010 V.Ivanchenko added member theGroundStateMass and inline
41//            method which allowing to compute this value once and use
42//            many times
[1340]43// 26.09.2010 V.Ivanchenko added number of protons, neutrons, proton holes
44//            and neutron holes as members of the class and Get/Set methods;
45//            removed not needed 'const'; removed old debug staff and unused
46//            private methods; add comments and reorder methods for
47//            better reading
[819]48
49#ifndef G4Fragment_h
50#define G4Fragment_h 1
51
52#include <vector>
53
54#include "globals.hh"
55#include "G4LorentzVector.hh"
56#include "G4ThreeVector.hh"
57#include "G4NucleiProperties.hh"
[1340]58//#include "G4ParticleTable.hh"
59//#include "G4IonTable.hh"
[819]60#include "Randomize.hh"
61#include "G4Proton.hh"
62#include "G4Neutron.hh"
63#include "G4HadTmpUtil.hh"
64
65class G4ParticleDefinition;
66
[1340]67class G4Fragment;     
[819]68typedef std::vector<G4Fragment*> G4FragmentVector;
69
70class G4Fragment 
71{
72public:
73
74  // ============= CONSTRUCTORS ==================
75
[1315]76  // Default constructor - obsolete
[819]77  G4Fragment();
78
79  // Destructor
80  ~G4Fragment();
81
82  // Copy constructor
83  G4Fragment(const G4Fragment &right);
84
[1315]85  // A,Z and 4-momentum - main constructor for fragment
[1340]86  G4Fragment(G4int A, G4int Z, const G4LorentzVector& aMomentum);
[819]87
[1340]88  // 4-momentum and pointer to G4particleDefinition (for gammas, e-)
89  G4Fragment(const G4LorentzVector& aMomentum, 
90             G4ParticleDefinition* aParticleDefinition);
[819]91
92  // ============= OPERATORS ==================
93   
94  const G4Fragment & operator=(const G4Fragment &right);
95  G4bool operator==(const G4Fragment &right) const;
96  G4bool operator!=(const G4Fragment &right) const;
97
98  friend std::ostream& operator<<(std::ostream&, const G4Fragment*);
99  friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
100
[1340]101  // ============= GENERAL METHODS ==================
[819]102
[1315]103  inline G4int GetZ_asInt() const;
104  inline G4int GetA_asInt() const;
105  inline void SetZandA_asInt(G4int Znew, G4int Anew);
[819]106 
[1315]107  inline G4double GetExcitationEnergy() const;
[1340]108
109  inline G4double GetGroundStateMass() const;
110   
111  inline G4double GetBindingEnergy() const;
[819]112 
[1315]113  inline const G4LorentzVector& GetMomentum() const;
114  inline void SetMomentum(const G4LorentzVector& value);
[819]115 
[1315]116  inline const G4ThreeVector& GetAngularMomentum() const;
117  inline void SetAngularMomentum(const G4ThreeVector& value);
[1340]118
119  // computation of mass for any Z and A
120  inline G4double ComputeGroundStateMass(G4int Z, G4int A) const;
121
122  // obsolete methods
123  inline G4double GetZ() const;
124  inline G4double GetA() const;
125  inline void SetZ(G4double value);
126  inline void SetA(G4double value);
[819]127 
[1340]128  // ============= METHODS FOR PRE-COMPOUND MODEL ===============
129
[1315]130  inline G4int GetNumberOfExcitons() const;
[819]131 
[1340]132  inline G4int GetNumberOfParticles() const;
133  inline G4int GetNumberOfCharged() const;
134  inline void SetNumberOfExcitedParticle(G4int valueTot, G4int valueP);
135
[1315]136  inline G4int GetNumberOfHoles() const;
[1340]137  inline G4int GetNumberOfChargedHoles() const;
138  inline void SetNumberOfHoles(G4int valueTot, G4int valueP=0);
[819]139 
[1340]140  // these methods will be removed in future
141  inline void SetNumberOfParticles(G4int value);
142  inline void SetNumberOfCharged(G4int value);
[819]143
[1340]144  // ============= METHODS FOR PHOTON EVAPORATION ===============
[819]145
[1340]146  inline G4int GetNumberOfElectrons() const;
147  inline void SetNumberOfElectrons(G4int value);
148
[1315]149  inline G4ParticleDefinition * GetParticleDefinition() const;
[1340]150  inline void SetParticleDefinition(G4ParticleDefinition * p);
[819]151
[1315]152  inline G4double GetCreationTime() const;
[1340]153  inline void SetCreationTime(G4double time);
[819]154
[1340]155  // ============= PRIVATE METHODS ==============================
[819]156
[1340]157private:
[819]158
[1340]159  void ExcitationEnergyWarning();
[819]160
[1340]161  void NumberOfExitationWarning(const G4String&);
[819]162
[1340]163  inline void CalculateExcitationEnergy();
[819]164
[1315]165  inline void CalculateGroundStateMass();
166
[819]167  // ============= DATA MEMBERS ==================
168
[1315]169  static G4int errCount;
170
171  G4int theA;
[819]172 
[1315]173  G4int theZ;
[819]174 
175  G4double theExcitationEnergy;
176
[1315]177  G4double theGroundStateMass;
178
[819]179  G4LorentzVector theMomentum;
180 
181  G4ThreeVector theAngularMomentum;
[1340]182
183  // Exciton model data members
[819]184 
185  G4int numberOfParticles;
186 
[1340]187  G4int numberOfCharged;
188
[819]189  G4int numberOfHoles;
190 
[1340]191  G4int numberOfChargedHoles;
[819]192
[1340]193  // Gamma evaporation data members
[819]194
[1340]195  G4int numberOfShellElectrons;
196
[819]197  G4ParticleDefinition * theParticleDefinition;
198 
199  G4double theCreationTime;
200
201};
202
[1340]203// ============= INLINE METHOD IMPLEMENTATIONS ===================
[819]204
[1340]205inline void G4Fragment::CalculateExcitationEnergy()
[819]206{
[1340]207  theExcitationEnergy = theMomentum.mag() - theGroundStateMass;
208  if(theExcitationEnergy < 0.0) { ExcitationEnergyWarning(); }
[819]209}
[1340]210         
211inline void G4Fragment::CalculateGroundStateMass() 
[819]212{
[1340]213  theGroundStateMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
[819]214}
215
[1315]216inline G4int G4Fragment::GetA_asInt() const
217{
218  return theA;
219}
220
221inline G4int G4Fragment::GetZ_asInt()  const
222{
223  return theZ;
224}
225
226inline void G4Fragment::SetZandA_asInt(G4int Znew, G4int Anew)
227{
228  theZ = Znew;
229  theA = Anew;
230  CalculateGroundStateMass();
231}
232
[819]233inline G4double G4Fragment::GetExcitationEnergy()  const
234{
235  return theExcitationEnergy;
236}
237
[1340]238inline G4double G4Fragment::GetGroundStateMass() const
239{
240  return theGroundStateMass; 
241}
242
243inline G4double G4Fragment::GetBindingEnergy() const
244{
245  return (theA-theZ)*CLHEP::neutron_mass_c2 + theZ*CLHEP::proton_mass_c2
246    - theGroundStateMass;
247}
248
[1315]249inline const G4LorentzVector& G4Fragment::GetMomentum()  const
[819]250{
251  return theMomentum;
252}
253
[1340]254inline void G4Fragment::SetMomentum(const G4LorentzVector& value)
255{
256  theMomentum = value;
257  CalculateExcitationEnergy();
258}
259
[1315]260inline const G4ThreeVector& G4Fragment::GetAngularMomentum()  const
[819]261{
262  return theAngularMomentum;
263}
264
[1315]265inline void G4Fragment::SetAngularMomentum(const G4ThreeVector& value)
[819]266{
267  theAngularMomentum = value;
268}
269
[1340]270inline G4double
271G4Fragment::ComputeGroundStateMass(G4int Z, G4int A) const
[819]272{
[1340]273  return G4NucleiProperties::GetNuclearMass(A, Z); 
[819]274}
275
[1340]276inline G4double G4Fragment::GetZ()  const
[819]277{
[1340]278  return G4double(theZ);
[819]279}
280
[1340]281inline G4double G4Fragment::GetA() const
[819]282{
[1340]283  return G4double(theA);
[819]284}
285
[1340]286inline void G4Fragment::SetZ(const G4double value)
[819]287{
[1340]288  theZ = G4lrint(value);
289  CalculateGroundStateMass();
[819]290}
291
[1340]292inline void G4Fragment::SetA(const G4double value)
293{
294  theA = G4lrint(value);
295  CalculateGroundStateMass();
296}
297
298inline G4int G4Fragment::GetNumberOfExcitons()  const
299{
300  return numberOfParticles + numberOfHoles;
301}
302
303inline G4int G4Fragment::GetNumberOfParticles()  const
304{
305  return numberOfParticles;
306}
307
[819]308inline G4int G4Fragment::GetNumberOfCharged()  const
309{
310  return numberOfCharged;
311}
312
[1340]313inline 
314void G4Fragment::SetNumberOfExcitedParticle(G4int valueTot, G4int valueP)
[819]315{
[1340]316  numberOfParticles = valueTot;
317  numberOfCharged = valueP;
318  if(valueTot < valueP)  { 
319    NumberOfExitationWarning("SetNumberOfExcitedParticle"); 
[819]320  }
321}
322
[1340]323inline G4int G4Fragment::GetNumberOfHoles()  const
[819]324{
[1340]325  return numberOfHoles;
[819]326}
327
[1340]328inline G4int G4Fragment::GetNumberOfChargedHoles()  const
[819]329{
[1340]330  return numberOfChargedHoles;
[819]331}
332
[1340]333inline void G4Fragment::SetNumberOfHoles(G4int valueTot, G4int valueP)
[819]334{
[1340]335  numberOfHoles = valueTot;
336  numberOfChargedHoles = valueP;
337  if(valueTot < valueP)  { 
338    NumberOfExitationWarning("SetNumberOfHoles"); 
339  }
[819]340}
341
[1340]342inline void G4Fragment::SetNumberOfParticles(G4int value)
[819]343{
[1340]344  numberOfParticles = value;
[819]345}
346
[1340]347inline void G4Fragment::SetNumberOfCharged(G4int value)
[819]348{
[1340]349  numberOfCharged = value;
350  if(value > numberOfParticles)  { 
351    NumberOfExitationWarning("SetNumberOfCharged"); 
352  }
[819]353}
354
[1340]355inline G4int G4Fragment::GetNumberOfElectrons() const
[819]356{
[1340]357  return numberOfShellElectrons;
[819]358}
[1315]359
[1340]360inline void G4Fragment::SetNumberOfElectrons(G4int value)
[1315]361{
[1340]362  numberOfShellElectrons = value;
[1315]363}
364
[1340]365inline 
366G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const
[1315]367{
[1340]368  return theParticleDefinition;
[1315]369}
[1340]370
371inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * p)
[819]372{
[1340]373  theParticleDefinition = p;
[819]374}
375
[1340]376inline G4double G4Fragment::GetCreationTime() const 
[1315]377{
[1340]378  return theCreationTime;
[1315]379}
[819]380
[1340]381inline void G4Fragment::SetCreationTime(G4double time)
382{
383  theCreationTime = time;
384}
385
[819]386#endif
387
388
Note: See TracBrowser for help on using the repository browser.