Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4QCandidate.hh
Go to the documentation of this file.
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$
28 //
29 // ---------------- G4QCandidate ----------------
30 // by Mikhail Kossov, Sept 1999.
31 // class header for Quasmon initiated Candidates used by the CHIPS Model
32 // ----------------------------------------------------------------------
33 // Short description: A candidate for hadronization. The candidates
34 // (hadrons or nuclear fragments) are competative, each quark of a
35 // Quasmon select which candidate to use for hadronization
36 // ------------------------------------------------------------------
37 
38 #ifndef G4QCandidate_h
39 #define G4QCandidate_h 1
40 
41 #include "G4QHadron.hh"
43 #include <algorithm>
44 
45 class G4QCandidate : public G4QHadron
46 {
47 public:
48  G4QCandidate(); // Default Constructor
49  G4QCandidate(G4int PDGcode); // Constructor by PDG Code
50  G4QCandidate(const G4QCandidate& right); // Copy Constructor by value
51  G4QCandidate(G4QCandidate* right); // Copy Constructor by pointer
52  ~G4QCandidate(); // Public Destructor
53  // Overloaded Operators
54  const G4QCandidate& operator=(const G4QCandidate& right);
55  G4bool operator==(const G4QCandidate &right) const;
56  G4bool operator!=(const G4QCandidate &right) const;
57  // Specific Selectors
58  G4QParentCluster* TakeParClust(G4int nPC);// Get pointer to ParentClust from ParClastVect
59  G4int GetPClustEntries() const; // Get a#of Parent Clusters in ParClastVector
60  G4bool GetPossibility() const; // Get possibility(true)/forbid(false) hadr/fr
61  G4bool GetParPossibility() const; // Get possibility(true)/forbidi(false) parent
62  G4double GetKMin() const; // Get k-minimal for the candidate
63  G4double GetEBMass() const; // Get bound mass in respect to Environment
64  G4double GetNBMass() const; // Get bound mass in respect to TotalNucleus
65  G4double GetDenseProbability() const; // Get dense-probability for the candidate
66  G4double GetPreProbability() const; // Get pre-probability for the candidate
67  G4double GetRelProbability() const; // Get the relative probility of hadronization
68  G4double GetIntegProbability() const; // Get integrated probability for randomization
69  G4double GetSecondRelProb() const; // Get 2nd relative probility of hadronization
70  G4double GetSecondIntProb() const; // Get 2nd integ. probability for randomization
71  // Specific Modifiers
72  void ClearParClustVector(); // Clear theParentClasterVector of theCandidate
73  void FillPClustVec(G4QParentCluster* pCl);// Set pointer to ParentClust in ParClastVector
74  void SetPossibility(G4bool choice); // Set possibility(true)/forbid(false) hadr/fr
75  void SetParPossibility(G4bool choice); // Set possibility(true)/forbid(false) parent
76  void SetKMin(G4double kmin); // Set k-minimal for the candidate
77  void SetDenseProbability(G4double prep); // Set dense-probability for the candidate
78  void SetPreProbability(G4double prep); // Set pre-probability for the candidate
79  void SetRelProbability(G4double relP); // Set the relative probility of hadronization
80  void SetIntegProbability(G4double intP); // Set integrated probability for randomization
81  void SetSecondRelProb(G4double relP); // Set 2nd relative probility of hadronization
82  void SetSecondIntProb(G4double intP); // Set 2nd integrProbability for randomization
83  void SetEBMass(G4double newMass); // Set mass bounded to Environment
84  void SetNBMass(G4double newMass); // Set mass bounded to Total Nucleus
85 
86 // Body
87 private:
88  G4bool possible; // permission/forbiding preFlag to be a hadron/fragment
89  G4bool parPossible; // permission/forbiding preFlag to be a parent
90  G4double kMin; // mu^2/2M (Q-case), ~BindingEnergy (Clust.-case)
91  G4double denseProbability; // a#of clusters of the type in dense region
92  G4double preProbability; // a#of clusters of the type or Q-suppression
93  G4double relativeProbability; // relative probability of hadronization
94  G4double integralProbability; // integrated probability of randomization
95  G4double secondRelProbability; // spare relative probability of hadronization
96  G4double secondIntProbability; // spare integrated probability of randomization
97  G4QParentClusterVector thePClusters; // vector of parent clusters for candid.-fragment
98  G4double EBMass; // Bound Mass of the cluster in the Environment
99  G4double NBMass; // Bound Mass of the cluster in the Total Nucleus
100 };
101 
102 inline G4bool G4QCandidate::operator==(const G4QCandidate &rhs) const {return this==&rhs;}
103 inline G4bool G4QCandidate::operator!=(const G4QCandidate &rhs) const {return this!=&rhs;}
104 
105 inline G4QParentCluster* G4QCandidate::TakeParClust(G4int nPC){return thePClusters[nPC];}
106 inline G4int G4QCandidate::GetPClustEntries() const {return thePClusters.size();}
107 inline G4bool G4QCandidate::GetPossibility() const {return possible;}
108 inline G4bool G4QCandidate::GetParPossibility() const {return parPossible;}
109 inline G4double G4QCandidate::GetKMin() const {return kMin;}
110 inline G4double G4QCandidate::GetEBMass() const {return EBMass;}
111 inline G4double G4QCandidate::GetNBMass() const {return NBMass;}
112 inline G4double G4QCandidate::GetDenseProbability() const {return denseProbability;}
113 inline G4double G4QCandidate::GetPreProbability() const {return preProbability;}
114 inline G4double G4QCandidate::GetRelProbability() const {return relativeProbability;}
115 inline G4double G4QCandidate::GetIntegProbability() const {return integralProbability;}
116 inline G4double G4QCandidate::GetSecondRelProb() const {return secondRelProbability;}
117 inline G4double G4QCandidate::GetSecondIntProb() const {return secondIntProbability;}
118 
120 {
121  std::for_each(thePClusters.begin(), thePClusters.end(), DeleteQParentCluster());
122  thePClusters.clear();
123 }
124 
126 {
127  thePClusters.push_back(pCl); // Fill new instance of PCl
128 }
131 inline void G4QCandidate::SetKMin(G4double kmin) {kMin=kmin;}
132 inline void G4QCandidate::SetDenseProbability(G4double prep) {denseProbability=prep;}
133 inline void G4QCandidate::SetPreProbability(G4double prep) {preProbability=prep;}
134 inline void G4QCandidate::SetRelProbability(G4double relP) {relativeProbability=relP;}
135 inline void G4QCandidate::SetIntegProbability(G4double intP) {integralProbability=intP;}
136 inline void G4QCandidate::SetSecondRelProb(G4double relP) {secondRelProbability=relP;}
137 inline void G4QCandidate::SetSecondIntProb(G4double intP) {secondIntProbability=intP;}
138 inline void G4QCandidate::SetEBMass(G4double newM) {EBMass=newM;}
139 inline void G4QCandidate::SetNBMass(G4double newM) {NBMass=newM;}
140 
141 #endif