00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // 00027 // $Id$ 00028 // 00029 // ---------------- G4QCandidate ---------------- 00030 // by Mikhail Kossov, Sept 1999. 00031 // class header for Quasmon initiated Candidates used by the CHIPS Model 00032 // ---------------------------------------------------------------------- 00033 // Short description: A candidate for hadronization. The candidates 00034 // (hadrons or nuclear fragments) are competative, each quark of a 00035 // Quasmon select which candidate to use for hadronization 00036 // ------------------------------------------------------------------ 00037 00038 #ifndef G4QCandidate_h 00039 #define G4QCandidate_h 1 00040 00041 #include "G4QHadron.hh" 00042 #include "G4QParentClusterVector.hh" 00043 #include <algorithm> 00044 00045 class G4QCandidate : public G4QHadron 00046 { 00047 public: 00048 G4QCandidate(); // Default Constructor 00049 G4QCandidate(G4int PDGcode); // Constructor by PDG Code 00050 G4QCandidate(const G4QCandidate& right); // Copy Constructor by value 00051 G4QCandidate(G4QCandidate* right); // Copy Constructor by pointer 00052 ~G4QCandidate(); // Public Destructor 00053 // Overloaded Operators 00054 const G4QCandidate& operator=(const G4QCandidate& right); 00055 G4bool operator==(const G4QCandidate &right) const; 00056 G4bool operator!=(const G4QCandidate &right) const; 00057 // Specific Selectors 00058 G4QParentCluster* TakeParClust(G4int nPC);// Get pointer to ParentClust from ParClastVect 00059 G4int GetPClustEntries() const; // Get a#of Parent Clusters in ParClastVector 00060 G4bool GetPossibility() const; // Get possibility(true)/forbid(false) hadr/fr 00061 G4bool GetParPossibility() const; // Get possibility(true)/forbidi(false) parent 00062 G4double GetKMin() const; // Get k-minimal for the candidate 00063 G4double GetEBMass() const; // Get bound mass in respect to Environment 00064 G4double GetNBMass() const; // Get bound mass in respect to TotalNucleus 00065 G4double GetDenseProbability() const; // Get dense-probability for the candidate 00066 G4double GetPreProbability() const; // Get pre-probability for the candidate 00067 G4double GetRelProbability() const; // Get the relative probility of hadronization 00068 G4double GetIntegProbability() const; // Get integrated probability for randomization 00069 G4double GetSecondRelProb() const; // Get 2nd relative probility of hadronization 00070 G4double GetSecondIntProb() const; // Get 2nd integ. probability for randomization 00071 // Specific Modifiers 00072 void ClearParClustVector(); // Clear theParentClasterVector of theCandidate 00073 void FillPClustVec(G4QParentCluster* pCl);// Set pointer to ParentClust in ParClastVector 00074 void SetPossibility(G4bool choice); // Set possibility(true)/forbid(false) hadr/fr 00075 void SetParPossibility(G4bool choice); // Set possibility(true)/forbid(false) parent 00076 void SetKMin(G4double kmin); // Set k-minimal for the candidate 00077 void SetDenseProbability(G4double prep); // Set dense-probability for the candidate 00078 void SetPreProbability(G4double prep); // Set pre-probability for the candidate 00079 void SetRelProbability(G4double relP); // Set the relative probility of hadronization 00080 void SetIntegProbability(G4double intP); // Set integrated probability for randomization 00081 void SetSecondRelProb(G4double relP); // Set 2nd relative probility of hadronization 00082 void SetSecondIntProb(G4double intP); // Set 2nd integrProbability for randomization 00083 void SetEBMass(G4double newMass); // Set mass bounded to Environment 00084 void SetNBMass(G4double newMass); // Set mass bounded to Total Nucleus 00085 00086 // Body 00087 private: 00088 G4bool possible; // permission/forbiding preFlag to be a hadron/fragment 00089 G4bool parPossible; // permission/forbiding preFlag to be a parent 00090 G4double kMin; // mu^2/2M (Q-case), ~BindingEnergy (Clust.-case) 00091 G4double denseProbability; // a#of clusters of the type in dense region 00092 G4double preProbability; // a#of clusters of the type or Q-suppression 00093 G4double relativeProbability; // relative probability of hadronization 00094 G4double integralProbability; // integrated probability of randomization 00095 G4double secondRelProbability; // spare relative probability of hadronization 00096 G4double secondIntProbability; // spare integrated probability of randomization 00097 G4QParentClusterVector thePClusters; // vector of parent clusters for candid.-fragment 00098 G4double EBMass; // Bound Mass of the cluster in the Environment 00099 G4double NBMass; // Bound Mass of the cluster in the Total Nucleus 00100 }; 00101 00102 inline G4bool G4QCandidate::operator==(const G4QCandidate &rhs) const {return this==&rhs;} 00103 inline G4bool G4QCandidate::operator!=(const G4QCandidate &rhs) const {return this!=&rhs;} 00104 00105 inline G4QParentCluster* G4QCandidate::TakeParClust(G4int nPC){return thePClusters[nPC];} 00106 inline G4int G4QCandidate::GetPClustEntries() const {return thePClusters.size();} 00107 inline G4bool G4QCandidate::GetPossibility() const {return possible;} 00108 inline G4bool G4QCandidate::GetParPossibility() const {return parPossible;} 00109 inline G4double G4QCandidate::GetKMin() const {return kMin;} 00110 inline G4double G4QCandidate::GetEBMass() const {return EBMass;} 00111 inline G4double G4QCandidate::GetNBMass() const {return NBMass;} 00112 inline G4double G4QCandidate::GetDenseProbability() const {return denseProbability;} 00113 inline G4double G4QCandidate::GetPreProbability() const {return preProbability;} 00114 inline G4double G4QCandidate::GetRelProbability() const {return relativeProbability;} 00115 inline G4double G4QCandidate::GetIntegProbability() const {return integralProbability;} 00116 inline G4double G4QCandidate::GetSecondRelProb() const {return secondRelProbability;} 00117 inline G4double G4QCandidate::GetSecondIntProb() const {return secondIntProbability;} 00118 00119 inline void G4QCandidate::ClearParClustVector() 00120 { 00121 std::for_each(thePClusters.begin(), thePClusters.end(), DeleteQParentCluster()); 00122 thePClusters.clear(); 00123 } 00124 00125 inline void G4QCandidate::FillPClustVec(G4QParentCluster* pCl) 00126 { 00127 thePClusters.push_back(pCl); // Fill new instance of PCl 00128 } 00129 inline void G4QCandidate::SetPossibility(G4bool choice) {possible=choice;} 00130 inline void G4QCandidate::SetParPossibility(G4bool choice) {parPossible=choice;} 00131 inline void G4QCandidate::SetKMin(G4double kmin) {kMin=kmin;} 00132 inline void G4QCandidate::SetDenseProbability(G4double prep) {denseProbability=prep;} 00133 inline void G4QCandidate::SetPreProbability(G4double prep) {preProbability=prep;} 00134 inline void G4QCandidate::SetRelProbability(G4double relP) {relativeProbability=relP;} 00135 inline void G4QCandidate::SetIntegProbability(G4double intP) {integralProbability=intP;} 00136 inline void G4QCandidate::SetSecondRelProb(G4double relP) {secondRelProbability=relP;} 00137 inline void G4QCandidate::SetSecondIntProb(G4double intP) {secondIntProbability=intP;} 00138 inline void G4QCandidate::SetEBMass(G4double newM) {EBMass=newM;} 00139 inline void G4QCandidate::SetNBMass(G4double newM) {NBMass=newM;} 00140 00141 #endif