// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // // ------------------------------------------------------------ // GEANT 4 class implementation file // // ---------------- G4FragmentingString ---------------- // by Gunter Folger, September 2001. // class for an excited string used in Fragmention // ------------------------------------------------------------ // G4FragmentingString #include "G4FragmentingString.hh" #include "G4ExcitedString.hh" //--------------------------------------------------------------------------------- //--------------------------------------------------------------------------------- G4FragmentingString::G4FragmentingString(const G4FragmentingString &old) { LeftParton=old.LeftParton; RightParton=old.RightParton; Ptleft=old.Ptleft; Ptright=old.Ptright; Pplus=old.Pplus; Pminus=old.Pminus; decaying=old.decaying; } //--------------------------------------------------------------------------------- G4FragmentingString::G4FragmentingString(const G4ExcitedString &excited) { LeftParton=excited.GetLeftParton()->GetDefinition(); RightParton=excited.GetRightParton()->GetDefinition(); Ptleft=excited.GetLeftParton()->Get4Momentum().vect(); Ptleft.setZ(0.); Ptright=excited.GetRightParton()->Get4Momentum().vect(); Ptright.setZ(0.); G4LorentzVector P=excited.Get4Momentum(); Pplus =P.e() + P.pz(); Pminus=P.e() - P.pz(); decaying=None; } //--------------------------------------------------------------------------------- G4FragmentingString::G4FragmentingString(const G4FragmentingString &old, G4ParticleDefinition * newdecay, const G4LorentzVector *momentum) { decaying=None; if ( old.decaying == Left ) { //G4cout<<" Left "<vect(); Ptleft.setZ(0.); //G4cout<<"Pt right "<vect(); Ptright.setZ(0.); LeftParton = old.LeftParton; Ptleft = old.Ptleft; //G4cout<<"Pt right "<e() + momentum->pz()); Pminus = old.Pminus - (momentum->e() - momentum->pz()); //G4double Eold=0.5 * (old.Pplus + old.Pminus); //G4double Enew=0.5 * (Pplus + Pminus); } //--------------------------------------------------------------------------------- G4FragmentingString::G4FragmentingString(const G4FragmentingString &old, // Uzhi G4ParticleDefinition * newdecay) // Uzhi { // Uzhi decaying=None; // Uzhi if ( old.decaying == Left ) // Uzhi { // Uzhi RightParton= old.RightParton; // Uzhi LeftParton = newdecay; // Uzhi } else if ( old.decaying == Right ) // Uzhi { // Uzhi RightParton = newdecay; // Uzhi LeftParton = old.LeftParton; // Uzhi } else // Uzhi { throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::G4FragmentingString: no decay Direction defined"); } } //--------------------------------------------------------------------------------- G4FragmentingString::~G4FragmentingString() {} //--------------------------------------------------------------------------------- void G4FragmentingString::SetLeftPartonStable() { theStableParton=GetLeftParton(); theDecayParton=GetRightParton(); decaying=Right; } //--------------------------------------------------------------------------------- void G4FragmentingString::SetRightPartonStable() { theStableParton=GetRightParton(); theDecayParton=GetLeftParton(); decaying=Left; } //--------------------------------------------------------------------------------- G4int G4FragmentingString::GetDecayDirection() const { if (decaying == Left ) return +1; else if (decaying == Right) return -1; else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::GetDecayDirection: decay side UNdefined!"); return 0; } //--------------------------------------------------------------------------------- G4bool G4FragmentingString::FourQuarkString() const { return LeftParton->GetParticleSubType()== "di_quark" && RightParton->GetParticleSubType()== "di_quark"; } //--------------------------------------------------------------------------------- G4bool G4FragmentingString::DecayIsQuark() { return theDecayParton->GetParticleSubType()== "quark"; } G4bool G4FragmentingString::StableIsQuark() { return theStableParton->GetParticleSubType()== "quark"; } //--------------------------------------------------------------------------------- G4ThreeVector G4FragmentingString::StablePt() { if (decaying == Left ) return Ptright; else if (decaying == Right ) return Ptleft; else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!"); return G4ThreeVector(); } G4ThreeVector G4FragmentingString::DecayPt() { if (decaying == Left ) return Ptleft; else if (decaying == Right ) return Ptright; else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!"); return G4ThreeVector(); } //--------------------------------------------------------------------------------- G4double G4FragmentingString::LightConePlus() { return Pplus; } G4double G4FragmentingString::LightConeMinus() { return Pminus; } G4double G4FragmentingString::LightConeDecay() { if (decaying == Left ) return Pplus; else if (decaying == Right ) return Pminus; else throw G4HadronicException(__FILE__, __LINE__, "G4FragmentingString::DecayPt: decay side UNdefined!"); return 0; } //--------------------------------------------------------------------------------- G4LorentzVector G4FragmentingString::Get4Momentum() const { G4LorentzVector momentum(Ptleft+Ptright,0); momentum.setPz(0.5*(Pplus-Pminus)); momentum.setE(0.5*(Pplus+Pminus)); return momentum; } G4double G4FragmentingString::Mass2() const { return Pplus*Pminus - (Ptleft+Ptright).mag2(); } G4double G4FragmentingString::Mass() const { return std::sqrt(this->Mass2()); } G4double G4FragmentingString::MassT2() const { return Pplus*Pminus; }