Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4LorentzConvertor.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 // $Id: G4LorentzConvertor.hh 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // 20100108 Michael Kelsey -- Use G4LorentzVector internally
29 // 20100120 M. Kelsey -- BUG FIX: scm_momentum should be G4ThreeVector
30 // 20100126 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
31 // 20100519 M. Kelsey -- Add interfaces to pass G4InuclParticles directly
32 // 20100616 M. Kelsey -- Report bullet and target four-momenta when set
33 // 20100915 M. Kelsey -- Move constructors to .cc file, add initializers
34 // 20110602 M. Kelsey -- Drop some unnecessary kinematics intermediates
35 
36 #ifndef G4LORENTZ_CONVERTOR_HH
37 #define G4LORENTZ_CONVERTOR_HH
38 
39 #include "globals.hh"
40 #include "G4LorentzVector.hh"
41 #include "G4ThreeVector.hh"
42 
43 class G4InuclParticle;
44 
46 public:
48 
49  G4LorentzConvertor(const G4LorentzVector& bmom, G4double bmass,
50  const G4LorentzVector& tmom, G4double tmass);
51 
52  G4LorentzConvertor(const G4InuclParticle* bullet,
53  const G4InuclParticle* target);
54 
55  void setVerbose(G4int vb=0) { verboseLevel = vb; }
56 
57  void setBullet(const G4InuclParticle* bullet);
58  void setTarget(const G4InuclParticle* target);
59 
60  void setBullet(const G4InuclParticle& bullet) { setBullet(&bullet); }
61  void setTarget(const G4InuclParticle& target) { setTarget(&target); }
62 
63  // Use correct four-vectors as input
64  void setBullet(const G4LorentzVector& bmom) {
65  bullet_mom = bmom;
66  if (verboseLevel > 3) printBullet();
67  }
68 
69  void setTarget(const G4LorentzVector& bmom) {
70  target_mom = bmom;
71  if (verboseLevel > 3) printTarget();
72  }
73 
74  // These functions "repair" input 4-vectors using specified mass
75  void setBullet(const G4LorentzVector& bmom, G4double bmass) {
76  bullet_mom.setVectM(bmom.vect(), bmass);
77  if (verboseLevel > 3) printBullet();
78  }
79 
80  void setTarget(const G4LorentzVector& tmom, G4double tmass) {
81  target_mom.setVectM(tmom.vect(), tmass);
82  if (verboseLevel > 3) printTarget();
83  }
84 
85  // Select reference frame for boosts, rotations, etc.
86  void toTheCenterOfMass();
87  void toTheTargetRestFrame();
88  void fillKinematics(); // Common calculations after either of above
89 
91 
92  // Four-vectors of bullet and target in last chosen reference frame
93  const G4LorentzVector& getBullet() const { return bullet_mom; }
94  const G4LorentzVector& getTarget() const { return target_mom; }
95 
97  G4double getTotalSCMEnergy() const { return ecm_tot; }
98  G4double getSCMMomentum() const { return scm_momentum.rho(); }
99  G4double getTRSMomentum() const;
100 
101  G4LorentzVector rotate(const G4LorentzVector& mom) const;
102 
104  const G4LorentzVector& mom) const;
105 
106  G4bool reflectionNeeded() const;
107 
108  G4bool trivial() const { return degenerated; }
109 
110  // Reporting functions for diagnostics
111  void printBullet() const;
112  void printTarget() const;
113 
114 private:
115  static const G4double small;
116 
120 
121  G4LorentzVector scm_momentum; // CM momentum relative to target/bullet
122  G4ThreeVector scm_direction; // Unit vector to reduce repeated calcs
123 
124  // Buffer variables for doing ::rotate() calculations
130 };
131 
132 #endif // G4LORENTZ_CONVERTOR_HH
const XML_Char * target
Definition: expat.h:268
void setTarget(const G4LorentzVector &bmom)
void setTarget(const G4LorentzVector &tmom, G4double tmass)
G4ThreeVector scm_direction
void setTarget(const G4InuclParticle &target)
void setBullet(const G4InuclParticle *bullet)
void setBullet(const G4InuclParticle &bullet)
void setVerbose(G4int vb=0)
G4double getTRSMomentum() const
G4LorentzVector rotate(const G4LorentzVector &mom) const
G4LorentzVector scm_momentum
G4LorentzVector bullet_mom
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static const G4double small
double rho() const
void setBullet(const G4LorentzVector &bmom)
const G4LorentzVector & getTarget() const
void setTarget(const G4InuclParticle *target)
G4bool reflectionNeeded() const
int G4int
Definition: G4Types.hh:78
G4double getSCMMomentum() const
G4bool trivial() const
const G4LorentzVector & getBullet() const
G4double getKinEnergyInTheTRS() const
Hep3Vector vect() const
G4double getTotalSCMEnergy() const
void setBullet(const G4LorentzVector &bmom, G4double bmass)
G4LorentzVector target_mom
void setVectM(const Hep3Vector &spatial, double mass)