source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/src/G4QCandidate.cc @ 962

Last change on this file since 962 was 962, checked in by garnier, 15 years ago

update processes

File size: 6.3 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//
26//
27// $Id: G4QCandidate.cc,v 1.34 2006/11/27 10:44:53 mkossov Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//      ---------------- G4QCandidate ----------------
31//             by Mikhail Kossov, Sept 1999.
32//      class for Quasmon initiated Candidates used by CHIPS Model
33// ------------------------------------------------------------------
34
35//#define debug
36
37#include "G4QCandidate.hh"
38#include <algorithm>
39
40G4QCandidate::G4QCandidate() : 
41  G4QHadron(),possible(false),parPossible(false),kMin(0),denseProbability(0.),
42  preProbability(0.),relativeProbability(0.),integralProbability(0.),
43  secondRelProbability(0.),secondIntProbability(0.),EBMass(0.),NBMass(0.)
44                           
45{
46}
47
48G4QCandidate::G4QCandidate(G4int PDGcode) :
49  G4QHadron(PDGcode),possible(false),parPossible(false),kMin(0),denseProbability(0.),
50  preProbability(0.),relativeProbability(0.),integralProbability(0.),
51  secondRelProbability(0.),secondIntProbability(0.),EBMass(0.),NBMass(0.)
52{
53#ifdef debug
54  G4cout<<"G4QCandidate::Constructor: PDG="<<PDGcode<<G4endl;
55#endif
56  G4LorentzVector cur4Mom(0.,0.,0.,0.);
57  G4QPDGCode QPDG(PDGcode);
58#ifdef debug
59  G4cout<<"G4QCandidate::Constructor: QPDG="<<QPDG<<G4endl;
60#endif
61  SetQPDG(QPDG);
62  G4double vacMass=QPDG.GetMass();
63#ifdef debug
64  G4cout<<"G4QCandidate::Constructor: M="<<vacMass<<G4endl;
65#endif
66  cur4Mom.setE(vacMass);
67  Set4Momentum(cur4Mom);
68  SetQC(QPDG.GetQuarkContent());
69}
70
71G4QCandidate::G4QCandidate(const G4QCandidate& right) : 
72  G4QHadron(&right)
73{
74  Set4Momentum         (right.Get4Momentum());
75  SetQPDG              (right.GetQPDG());
76  SetQC                (right.GetQC());
77  SetNFragments        (right.GetNFragments());
78  possible            = right.possible;
79  parPossible         = right.parPossible;
80  kMin                = right.kMin;
81  denseProbability    = right.denseProbability;
82  preProbability      = right.preProbability;
83  relativeProbability = right.relativeProbability;
84  integralProbability = right.integralProbability;
85  secondRelProbability= right.secondRelProbability;
86  secondIntProbability= right.secondIntProbability;
87  EBMass = right.EBMass;
88  NBMass = right.NBMass;
89  // thePClusters
90  G4int nParCl        = right.thePClusters.size();
91  if(nParCl) for(G4int ip=0; ip<nParCl; ip++)
92  {
93    G4QParentCluster* curPC = new G4QParentCluster(right.thePClusters[ip]);
94    thePClusters.push_back(curPC);
95  }
96}
97
98G4QCandidate::G4QCandidate(G4QCandidate* right)
99{
100  Set4Momentum         (right->Get4Momentum());
101  SetQPDG              (right->GetQPDG());
102  SetQC                (right->GetQC());
103  SetNFragments        (right->GetNFragments());
104  possible            = right->possible;
105  parPossible         = right->parPossible;
106  kMin                = right->kMin;
107  denseProbability    = right->denseProbability;
108  preProbability      = right->preProbability;
109  relativeProbability = right->relativeProbability;
110  integralProbability = right->integralProbability;
111  secondRelProbability= right->secondRelProbability;
112  secondIntProbability= right->secondIntProbability;
113  EBMass = right->EBMass;
114  NBMass = right->NBMass;
115  // thePClusters
116  G4int nParCl        = right->thePClusters.size();
117  if(nParCl) for(G4int ip=0; ip<nParCl; ip++)
118  {
119    G4QParentCluster* curPC = new G4QParentCluster(right->thePClusters[ip]);
120    thePClusters.push_back(curPC);
121  }
122}
123
124G4QCandidate::~G4QCandidate()
125{
126#ifdef debug
127  G4cout<<"~G4QCandidate: before thePClusters nC="<<thePClusters.size()<<G4endl;
128#endif
129  std::for_each(thePClusters.begin(), thePClusters.end(), DeleteQParentCluster());
130#ifdef debug
131  G4cout<<"~G4QCandidate: === DONE ==="<<G4endl;
132#endif
133}
134
135// Assignment operator
136const G4QCandidate& G4QCandidate::operator=(const G4QCandidate &right)
137{
138  if(this != &right)                          // Beware of self assignment
139  {
140    Set4Momentum         (right.Get4Momentum());
141    SetQPDG              (right.GetQPDG());
142    SetQC                (right.GetQC());
143    SetNFragments        (right.GetNFragments());
144    possible            = right.possible;
145    parPossible         = right.parPossible;
146    kMin                = right.kMin;
147    denseProbability    = right.denseProbability;
148    preProbability      = right.preProbability;
149    relativeProbability = right.relativeProbability;
150    integralProbability = right.integralProbability;
151    secondRelProbability= right.secondRelProbability;
152    secondIntProbability= right.secondIntProbability;
153    EBMass = right.EBMass;
154    NBMass = right.NBMass;
155    // thePClusters (Vector)
156    G4int nParCl        = right.thePClusters.size();
157    if(nParCl) for(G4int ip=0; ip<nParCl; ip++)
158    {
159      G4QParentCluster* curPC = new G4QParentCluster(right.thePClusters[ip]);
160      thePClusters.push_back(curPC);
161    }
162  }
163  return *this;
164}
Note: See TracBrowser for help on using the repository browser.