source: trunk/source/processes/hadronic/models/cascade/cascade/include/G4LorentzConvertor.hh

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

update ti head

File size: 4.8 KB
Line 
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// $Id: G4LorentzConvertor.hh,v 1.17 2010/09/15 20:16:16 mkelsey Exp $
26// Geant4 tag: $Name: hadr-casc-V09-03-85 $
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
35#ifndef G4LORENTZ_CONVERTOR_HH
36#define G4LORENTZ_CONVERTOR_HH
37
38#include "globals.hh"
39#include "G4LorentzVector.hh"
40#include "G4ThreeVector.hh"
41
42class G4InuclParticle;
43
44class G4LorentzConvertor {
45public:
46  G4LorentzConvertor();
47
48  G4LorentzConvertor(const G4LorentzVector& bmom, G4double bmass, 
49                     const G4LorentzVector& tmom, G4double tmass);
50
51  G4LorentzConvertor(const G4InuclParticle* bullet, 
52                     const G4InuclParticle* target);
53
54  void setVerbose(G4int vb=0) { verboseLevel = vb; }
55
56  void setBullet(const G4InuclParticle* bullet);
57  void setTarget(const G4InuclParticle* target);
58
59  void setBullet(const G4InuclParticle& bullet) { setBullet(&bullet); }
60  void setTarget(const G4InuclParticle& target) { setTarget(&target); }
61
62  // Use correct four-vectors as input
63  void setBullet(const G4LorentzVector& bmom) {
64    bullet_mom = bmom;
65    if (verboseLevel > 3) printBullet();
66  }
67
68  void setTarget(const G4LorentzVector& bmom) {
69    target_mom = bmom;
70    if (verboseLevel > 3) printTarget();
71  }
72
73  // These functions "repair" input 4-vectors using specified mass
74  void setBullet(const G4LorentzVector& bmom, G4double bmass) {
75    bullet_mom.setVectM(bmom.vect(), bmass);
76    if (verboseLevel > 3) printBullet();
77  }
78 
79  void setTarget(const G4LorentzVector& tmom, G4double tmass) {
80    target_mom.setVectM(tmom.vect(), tmass);
81    if (verboseLevel > 3) printTarget();
82  }
83
84  void toTheCenterOfMass();
85  void toTheTargetRestFrame(); 
86
87  G4LorentzVector backToTheLab(const G4LorentzVector& mom) const;
88
89  // Four-vectors of bullet and target in last chosen reference frame
90  const G4LorentzVector& getBullet() const { return bullet_mom; }
91  const G4LorentzVector& getTarget() const { return target_mom; }
92 
93  G4double getKinEnergyInTheTRS() const;
94  G4double getTotalSCMEnergy() const { return ecm_tot; }
95  G4double getSCMMomentum() const { return scm_momentum.rho(); }
96  G4double getTRSMomentum() const;
97
98  G4LorentzVector rotate(const G4LorentzVector& mom) const; 
99
100  G4LorentzVector rotate(const G4LorentzVector& mom1,
101                         const G4LorentzVector& mom) const; 
102
103  G4bool reflectionNeeded() const; 
104
105  G4bool trivial() const { return degenerated; }
106
107  // Reporting functions for diagnostics
108  void printBullet() const;
109  void printTarget() const;
110
111private: 
112  static const G4double small;
113
114  G4int verboseLevel;
115  G4LorentzVector bullet_mom;
116  G4LorentzVector target_mom;
117
118  G4LorentzVector scm_momentum;         // CM momentum relative to target/bullet
119
120  // Buffer variables for doing ::rotate() calculations
121  G4ThreeVector velocity;
122  G4double gamma;
123  G4double v2;
124  G4double ecm_tot;
125  G4double ga;
126  G4double gb;
127  G4double gbpp;
128  G4double gapp;
129  G4bool degenerated;
130};       
131
132#endif // G4LORENTZ_CONVERTOR_HH
Note: See TracBrowser for help on using the repository browser.