// // ******************************************************************** // * 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 (SoftHadrons) // ---------------------------------------------------------------- // //#define debug //#define trapdebug //#define pdebug //#define ppdebug #include "G4StringChipsParticleLevelInterface.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 #ifdef hdebug_SCPLI const G4int G4StringChipsParticleLevelInterface::nbh=200; G4double G4StringChipsParticleLevelInterface::bhmax=20.; G4double G4StringChipsParticleLevelInterface::ehmax=20.; G4double G4StringChipsParticleLevelInterface::bhdb=0.; G4double G4StringChipsParticleLevelInterface::ehde=0.; G4double G4StringChipsParticleLevelInterface::toth=0.; G4int G4StringChipsParticleLevelInterface::bover=0; G4int G4StringChipsParticleLevelInterface::eover=0; G4int* G4StringChipsParticleLevelInterface::bhis = new G4int[G4StringChipsParticleLevelInterface::nbh]; G4int* G4StringChipsParticleLevelInterface::ehis = new G4int[G4StringChipsParticleLevelInterface::nbh]; #endif G4StringChipsParticleLevelInterface::G4StringChipsParticleLevelInterface() { #ifdef debug G4cout<<"G4StringChipsParticleLevelInterface::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 the cut-off for calculating the nuclear radius in percent"<> theInnerCoreDensityCut; } } G4HadFinalState* G4StringChipsParticleLevelInterface:: ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus) { #ifdef debug G4cout<<"G4StringChipsParticleLevelInterface::ApplyYourself is called"<GetPDGMass(); static const G4double mNeut=G4Neutron::Neutron()->GetPDGMass(); static const G4double mLamb=G4Lambda::Lambda()->GetPDGMass(); static const G4double mKChg=G4KaonPlus::KaonPlus()->GetPDGMass(); static const G4double mKZer=G4KaonZero::KaonZero()->GetPDGMass(); static const G4double mPiCh=G4PionMinus::PionMinus()->GetPDGMass(); static const G4int pcl=4; // clusterization parameter for Energy Flow static const G4QContent ProtQC(1,2,0,0,0,0); static const G4QContent NeutQC(2,1,0,0,0,0); static const G4QContent LambQC(1,1,1,0,0,0); static const G4QContent KPlsQC(0,1,0,0,0,1); static const G4QContent KMinQC(0,0,1,0,1,0); static const G4QContent AKZrQC(1,0,0,0,0,1); static const G4QContent KZerQC(1,0,0,0,0,1); static const G4QContent PiMiQC(1,0,0,0,1,0); static const G4QContent PiPlQC(0,1,0,1,0,0); #ifdef debug G4cout<<"G4StringChipsParticleLevelInterface::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__, // "G4StringChipsParticleLevelInterface: 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); //G4double radius2 = theNucleus->GetNuclearRadius(theInnerCoreDensityCut*perCent); radius2 *= radius2; G4double pathlength = 0.; #ifdef ppdebug G4cout<<"G4StringChipsParticleLevelInterface::Propagate: r="<operator[](secondary)->GetDefinition()->GetPDGCharge()<<" " << theSecondaries->operator[](secondary)->GetDefinition()->GetPDGEncoding()<<" " << a4Mom <Get4Momentum().rapidity()<GetDefinition()->GetPDGCharge()<<" " << (*current).second->GetDefinition()->GetPDGEncoding()<<" " << (*current).second->Get4Momentum() <0 collection ---------- G4QContent EnFlowQC(0,0,0,0,0,0); G4LorentzVector EnFlow4M(0.,0.,0.,0.); //G4bool empty=true; G4int barys=0; G4int stras=0; G4int chars=0; for(G4int hp = theContents.size()-1; hp>=0; hp--) { G4QContent curQC=theContents[hp]; G4int baryn = curQC.GetBaryonNumber(); G4int stran = curQC.GetStrangeness(); G4int charg = curQC.GetCharge(); EnFlowQC += curQC; // Keep collecting energy flow EnFlow4M += *(theMomenta[hp]); barys += baryn; // Collected baryon number stras += stran; // Collected strangeness chars += charg; // Collected charge //empty = false; } if(barys>pcl) // Split in two or more parts (to survive!) { G4int nprt=(barys-1)/pcl+1; // Number of parts (pcl=4: 2:5-8,3:9-12...) G4int curb=barys; while (nprt>0) { nprt--; // One part is going to be created G4int brnm=pcl; // Baryon number of splitting part curb-=brnm; // The residual baryon number G4double prtM=0.; // The resulting GS mass of the part G4double resM=0.; // The resulting GS mass of the residual G4QContent prtQC(0,0,0,0,0,0); // The resulting Quark Content of the part G4int strm=0; // Max strangeness per part (stras=0) if(stras>0) strm=(stras-1)/nprt+1; // Max strangeness per part (stras>0) else if(stras<0) strm=(stras+1)/nprt-1; // Max strangeness per part (stras<0) G4int chgm=0; // Max charge per part (chars=0) if(stras>0) chgm=(chars-1)/nprt+1; // Max strangeness per part (chars>0) else if(stras<0) chgm=(chars+1)/nprt-1; // Max strangeness per part (chars<0) // ---> calculate proposed separated part //@@ Convert it to a CHIPS function (Which class? G4QH::Conctruct?) if(!strm) // --> The total strangness = 0 (n/p/pi-) { if(chgm<0) // (n/pi-) { prtM=(-chgm)*mPiCh+brnm*mNeut; prtQC=(-chgm)*PiMiQC+brnm*NeutQC; } else // (n/p) { prtM=chgm*mProt+(brnm-chgm)*mNeut; prtQC=chgm*ProtQC+(brnm-chgm)*NeutQC; } } else if(strm>=brnm) // ---> BigPositiveStrangeness(L/Pi+/K0/K-) { G4int stmb=strm-brnm; if(chgm<0) // (L/K-/K0) { prtM=(-chgm)*mKChg+brnm*mLamb+std::abs(stmb+chgm)*mKZer; prtQC=(-chgm)*KMinQC+brnm*LambQC; if(stmb>-chgm) prtQC+=(stmb+chgm)*KZerQC; else if(stmb<-chgm) prtQC+=(-stmb-chgm)*AKZrQC; } else // (L/K0/pi+) { prtM=chgm*mPiCh+(strm-brnm)*mKZer+brnm*mLamb; prtQC=chgm*PiPlQC+(strm-brnm)*KZerQC+brnm*LambQC; } } else if(strm>0) // ---> PositiveStrangeness=bmst) // (L/p/Pi+) { prtM=(chgm-bmst)*mPiCh+strm*mLamb+bmst*mProt; prtQC=(chgm-bmst)*PiPlQC+strm*LambQC+bmst*ProtQC; } else // ch NegativeStrangeness (N/K+/aK0/Pi-) { G4int bmst=brnm-strm; if(chgm>=bmst) // (K+/p/Pi+) { prtM=(-strm)*mKChg+brnm*mProt+(chgm-bmst)*mPiCh; prtQC=(-strm)*KPlsQC+brnm*ProtQC+(chgm-bmst)*PiPlQC; } else if(chgm>=-strm) // (K+/p/n) { prtM=(-strm)*mKChg+chgm*mProt+(brnm-chgm)*mNeut; prtQC=(-strm)*KPlsQC+chgm*ProtQC+(brnm-chgm)*NeutQC; } else if(chgm>=0) // (K+/aK0/n) { prtM=chgm*mKChg+(-chgm-strm)*mKZer+brnm*mNeut; prtQC=chgm*KPlsQC+(-chgm-strm)*AKZrQC+brnm*NeutQC; } else // ch<0 (aK0/n/Pi-) { prtM=(-strm)*mKChg+(-chgm)*mPiCh+brnm*mNeut; prtQC=(-strm)*KPlsQC+(-chgm)*PiMiQC+brnm*NeutQC; } } EnFlowQC-=prtQC; chgm=chars-chgm; // Just to keep the same notation strm=stras-strm; brnm=curb; if(!strm) // --> The total strangness = 0 (n/p/pi-) { if(chgm<0) resM=(-chgm)*mPiCh+brnm*mNeut; else resM=chgm*mProt+(brnm-chgm)*mNeut; } else if(strm>=brnm) // ---> BigPositiveStrangeness(L/Pi+/K0/K-) { G4int stmb=strm-brnm; if(chgm<0) resM=(-chgm)*mKChg+brnm*mLamb+std::abs(stmb+chgm)*mKZer; else resM=chgm*mPiCh+(strm-brnm)*mKZer+brnm*mLamb; } else if(strm>0) // ---> PositiveStrangeness=bmst) resM=(chgm-bmst)*mPiCh+strm*mLamb+bmst*mProt; else resM=chgm*mProt+strm*mLamb+(bmst-chgm)*mNeut; } else // ---> NegativeStrangeness (N/K+/aK0/Pi-) { G4int bmst=brnm-strm; if (chgm>=bmst) resM=(-strm)*mKChg+brnm*mProt+(chgm-bmst)*mPiCh; else if(chgm>=-strm) resM=(-strm)*mKChg+chgm*mProt+(brnm-chgm)*mNeut; else if(chgm>=0) resM=chgm*mKChg+(-chgm-strm)*mKZer+brnm*mNeut; else resM=(-strm)*mKChg+(-chgm)*mPiCh+brnm*mNeut; } G4LorentzVector prt4M=(prtM/(prtM+resM))*EnFlow4M; EnFlow4M-=prt4M; EnFlowQC-=prtQC; G4QHadron* aHadron = new G4QHadron(prtQC, prt4M); projHV.push_back(aHadron); if(nprt==1) { G4QHadron* fHadron = new G4QHadron(EnFlowQC, EnFlow4M); projHV.push_back(fHadron); nprt=0; } #ifdef debug G4cout<<"G4StringChipsParticleLevelInterface::Propagate: nprt="<0 collection ---------- too forward //G4QContent EnFlowQC(0,0,0,0,0,0); //G4LorentzVector EnFlow4M(0.,0.,0.,0.); //G4bool empty=true; //G4int barys=0; //for(G4int hp = theContents.size()-1; hp>=0; hp--) //{ // G4QContent curQC=theContents[hp]; // G4int baryn = curQC.GetBaryonNumber(); // if(baryn>0 && barys>0 && !empty) // New baryon and b-positive collection -> fill // { // G4QHadron* aHadron = new G4QHadron(EnFlowQC, EnFlow4M); // projHV.push_back(aHadron); // EnFlowQC = G4QContent(0,0,0,0,0,0); // EnFlow4M = G4LorentzVector(0.,0.,0.,0.); // barys=0; // empty=true; // } // EnFlowQC += curQC; // Keep collecting energy flow // EnFlow4M += *(theMomenta[hp]); // barys += baryn; // empty = false; //} //if(!empty) //{ // G4QHadron* aHadron = new G4QHadron(EnFlowQC, EnFlow4M); // projHV.push_back(aHadron); //} // End of Energy Flow One Quasmon for each ---------------------------- } else // Never come here (!?) { unsigned int hp; for(hp=0; hpGet4Momentum()<<" "<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()<