// // ******************************************************************** // * 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. * // ******************************************************************** // // Short description: Interface of QGSC to CHIPS (EnergyFlow) // ---------------------------------------------------------------------------- //#define debug //#define pdebug #include "G4QStringChipsParticleLevelInterface.hh" #include "globals.hh" #include #include #include #include "G4KineticTrackVector.hh" #include "G4Nucleon.hh" #include "G4Proton.hh" #include "G4Neutron.hh" #include "G4LorentzRotation.hh" #include "G4HadronicException.hh" // #define CHIPSdebug // #define CHIPSdebug_1 G4QStringChipsParticleLevelInterface::G4QStringChipsParticleLevelInterface() { #ifdef debug G4cout<<"G4QStringChipsParticleLevelInterface::Constructor is called"<> theEnergyLossPerFermi; theEnergyLossPerFermi *= GeV; G4cout << "Please enter nop"<> nop; G4cout << "Please enter the fractionOfSingleQuasiFreeNucleons"<> fractionOfSingleQuasiFreeNucleons; G4cout << "Please enter the fractionOfPairedQuasiFreeNucleons"<> fractionOfPairedQuasiFreeNucleons; G4cout << "Please enter the clusteringCoefficient"<> clusteringCoefficient; G4cout << "Please enter the temperature"<> temperature; G4cout << "Please enter halfTheStrangenessOfSee"<> halfTheStrangenessOfSee; G4cout << "Please enter the etaToEtaPrime"<> etaToEtaPrime; G4cout << "Please enter the fusionToExchange"<> fusionToExchange; G4cout << "Please enter cut-off for calculating the nuclear radius in percent"<> theInnerCoreDensityCut; } } G4HadFinalState* G4QStringChipsParticleLevelInterface:: ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus) { #ifdef debug G4cout<<"G4QStringChipsParticleLevelInterface::ApplyYourself is called"<GetPDGMass(); static const G4double mNeut=G4Neutron::Neutron()->GetPDGMass(); static const G4double mLamb=G4Lambda::Lambda()->GetPDGMass(); #ifdef debug G4cout<<"G4QStringChipsParticleLevelInterface::Propagate is called"<size() == 1) { G4ReactionProductVector* theFastResult = new G4ReactionProductVector; G4ReactionProduct* theFastSec; theFastSec = new G4ReactionProduct((*theSecondaries)[0]->GetDefinition()); G4LorentzVector current4Mom = (*theSecondaries)[0]->Get4Momentum(); theFastSec->SetTotalEnergy(current4Mom.t()); theFastSec->SetMomentum(current4Mom.vect()); theFastResult->push_back(theFastSec); return theFastResult; //throw G4HadronicException(__FILE__,__LINE__, // "G4QStringChipsParticleLevelInterface: Only one particle from String models!"); } // target properties needed in constructor of quasmon, and for boosting to // target rest frame // remove all nucleons already involved in STRING interaction, to make the ResidualTarget theNucleus->StartLoop(); G4Nucleon * aNucleon; G4int resA = 0; G4int resZ = 0; G4ThreeVector hitMomentum(0,0,0); G4double hitMass = 0; unsigned int hitCount = 0; while((aNucleon = theNucleus->GetNextNucleon())) { if(!aNucleon->AreYouHit()) { resA++; // Collect A of the ResidNuc resZ+=G4int (aNucleon->GetDefinition()->GetPDGCharge()); // Collect Z of the ResidNuc } else { hitMomentum += aNucleon->GetMomentum().vect(); // Sum 3-mom of StringHadr's hitMass += aNucleon->GetMomentum().m(); // Sum masses of StringHadrs hitCount ++; // Calculate STRING hadrons } } G4int targetPDGCode = 90000000 + 1000*resZ + (resA-resZ); // PDG of theResidualNucleus G4double targetMass=mNeut; if (!resZ) // Nucleus of only neutrons { if (resA>1) targetMass*=resA; } else targetMass=G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ) ->GetPDGMass(); G4double targetEnergy = std::sqrt(hitMomentum.mag2()+targetMass*targetMass); // !! @@ Target should be at rest: hitMomentum=(0,0,0) @@ !! M.K. (go to this system) G4LorentzVector targ4Mom(-1.*hitMomentum, targetEnergy); // Calculate the mean energy lost std::pair theImpact = theNucleus->RefetchImpactXandY(); G4double impactX = theImpact.first; G4double impactY = theImpact.second; G4double impactPar2 = impactX*impactX + impactY*impactY; G4double radius2 = theNucleus->GetNuclearRadius(theInnerCoreDensityCut*perCent); radius2 *= radius2; G4double pathlength = 0.; #ifdef pdebug G4cout<<"G4QStringChipsParticleLevelInterface::Propagate: r="<operator[](secondary)->GetDefinition()->GetPDGCharge()<<" " << theSecondaries->operator[](secondary)->GetDefinition()->GetPDGEncoding()<<" " << a4Mom <Get4Momentum().rapidity()<GetDefinition()->GetPDGCharge()<<" " << (*current).second->GetDefinition()->GetPDGEncoding()<<" " << (*current).second->Get4Momentum() <Get4Momentum()<<" "<GetPDGCode()<size()<GetPDGCharge()<<" " << theDefinition->GetPDGEncoding()<<" " << output->operator[](particle)->Get4Momentum() <operator[](particle); } else { if(resA>0) { G4ParticleDefinition* theDefinition = G4Neutron::Neutron(); if(resA==1) // The residual nucleus at rest must be added to conserve BaryN & Charge { if(resZ == 1) theDefinition = G4Proton::Proton(); } else theDefinition = G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ); theSec = new G4ReactionProduct(theDefinition); theSec->SetTotalEnergy(theDefinition->GetPDGMass()); theSec->SetMomentum(G4ThreeVector(0.,0.,0.)); theResult->push_back(theSec); if(!resZ && resA>0) for(G4int ni=1; niSetTotalEnergy(theDefinition->GetPDGMass()); theSec->SetMomentum(G4ThreeVector(0.,0.,0.)); theResult->push_back(theSec); } } } delete output; #ifdef CHIPSdebug G4cout << "Number of particles"<size()<