Geant4 Cross Reference |
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 // Historic fragment from M.Komogorov; clean-up still necessary @@@ 27 28 #include "G4ExcitedStringDecay.hh" 29 #include "G4SystemOfUnits.hh" 30 #include "G4KineticTrack.hh" 31 #include "G4LundStringFragmentation.hh" 32 #include "G4HadronicInteractionRegistry.hh" 33 34 #include "G4SampleResonance.hh" 35 36 //#define debug_G4ExcitedStringDecay 37 //#define debug_G4ExcitedStringCorr 38 39 G4ExcitedStringDecay::G4ExcitedStringDecay(G4VLongitudinalStringDecay* ptr) 40 : G4VStringFragmentation(), theStringDecay(ptr) 41 { 42 if(!ptr) { 43 G4HadronicInteraction* p = 44 G4HadronicInteractionRegistry::Instance()->FindModel("LundStringFragmentation"); 45 theStringDecay = static_cast<G4VLongitudinalStringDecay*>(p); 46 if(!theStringDecay) { theStringDecay = new G4LundStringFragmentation(); } 47 } 48 SetModelName(theStringDecay->GetModelName()); 49 } 50 51 G4ExcitedStringDecay::~G4ExcitedStringDecay() 52 {} 53 54 G4KineticTrackVector *G4ExcitedStringDecay::FragmentString(const G4ExcitedString &theString) 55 { 56 return theStringDecay->FragmentString(theString); 57 } 58 59 G4KineticTrackVector *G4ExcitedStringDecay::FragmentStrings(const G4ExcitedStringVector * theStrings) 60 { 61 G4LorentzVector KTsum(0.,0.,0.,0.); 62 63 #ifdef debug_G4ExcitedStringDecay 64 G4cout<<G4endl; 65 G4cout<<"--------------------------- G4ExcitedStringDecay ----------------------"<<G4endl; 66 G4cout<<"Hadronization of Excited Strings: theStrings->size() "<<theStrings->size()<<G4endl; 67 #endif 68 69 for ( unsigned int astring=0; astring < theStrings->size(); astring++) 70 // for ( unsigned int astring=0; astring < 1; astring++) 71 { 72 // G4cout<<"theStrings->operator[](astring)->IsExcited() "<<" "<<astring<<" "<<theStrings->operator[](astring)->IsExcited()<<G4endl; 73 if ( theStrings->operator[](astring)->IsExcited() ) 74 {KTsum+= theStrings->operator[](astring)->Get4Momentum();} 75 else {KTsum+=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum();} 76 } 77 78 G4LorentzRotation toCms( -1 * KTsum.boostVector() ); 79 G4LorentzRotation toLab(toCms.inverse()); 80 G4LorentzVector Ptmp; 81 KTsum=G4LorentzVector(0.,0.,0.,0.); 82 83 for ( unsigned int astring=0; astring < theStrings->size(); astring++) 84 // for ( unsigned int astring=0; astring < 1; astring++) 85 { 86 if ( theStrings->operator[](astring)->IsExcited() ) 87 { 88 Ptmp=toCms * theStrings->operator[](astring)->GetLeftParton()->Get4Momentum(); 89 theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp); 90 91 Ptmp=toCms * theStrings->operator[](astring)->GetRightParton()->Get4Momentum(); 92 theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp); 93 94 KTsum+= theStrings->operator[](astring)->Get4Momentum(); 95 } 96 else 97 { 98 Ptmp=toCms * theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum(); 99 theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp); 100 KTsum+= theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum(); 101 } 102 } 103 104 G4SampleResonance BrW; 105 const G4ParticleDefinition* TrackDefinition=0; 106 107 G4KineticTrackVector * theResult = new G4KineticTrackVector; 108 G4int attempts(0); 109 G4bool success=false; 110 G4bool NeedEnergyCorrector=false; 111 do { 112 #ifdef debug_G4ExcitedStringDecay 113 G4cout<<"New try No "<<attempts<<" to hadronize strings"<<G4endl; 114 #endif 115 116 std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack()); 117 theResult->clear(); 118 119 attempts++; 120 121 G4LorentzVector KTsecondaries(0.,0.,0.,0.); 122 NeedEnergyCorrector=false; 123 124 for ( unsigned int astring=0; astring < theStrings->size(); astring++) 125 // for ( unsigned int astring=0; astring < 1; astring++) // Proj. Last Str. Decay for FTF 126 // for ( unsigned int astring=1; astring < 2; astring++) // Proj. Last Str. Decay for QGS 127 // for ( unsigned int astring=0; astring < 1; astring++) // For testing purposes 128 { 129 #ifdef debug_G4ExcitedStringDecay 130 G4cout<<"String No "<<astring+1<<" Excited? "<<theStrings->operator[](astring)->IsExcited()<<G4endl; 131 132 G4cout<<"String No "<<astring+1<<" 4Momentum "<<theStrings->operator[](astring)->Get4Momentum() 133 <<" "<<theStrings->operator[](astring)->Get4Momentum().mag()<<G4endl; 134 #endif 135 136 G4KineticTrackVector * generatedKineticTracks = NULL; 137 if ( theStrings->operator[](astring)->IsExcited() ) 138 { 139 #ifdef debug_G4ExcitedStringDecay 140 G4cout<<"Fragment String with partons: " 141 <<theStrings->operator[](astring)->GetLeftParton()->GetPDGcode() <<" " 142 <<theStrings->operator[](astring)->GetRightParton()->GetPDGcode()<<" " 143 <<"Direction "<<theStrings->operator[](astring)->GetDirection()<<G4endl; 144 #endif 145 generatedKineticTracks=FragmentString(*theStrings->operator[](astring)); 146 #ifdef debug_G4ExcitedStringDecay 147 G4cout<<"(G4ExcitedStringDecay) Number of produced hadrons = " 148 <<generatedKineticTracks->size()<<G4endl; 149 #endif 150 } else { 151 #ifdef debug_G4ExcitedStringDecay 152 G4cout<<" GetTrack from the String"<<G4endl; 153 #endif 154 G4LorentzVector Mom=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum(); 155 G4KineticTrack * aTrack= new G4KineticTrack( 156 theStrings->operator[](astring)->GetKineticTrack()->GetDefinition(), 157 theStrings->operator[](astring)->GetKineticTrack()->GetFormationTime(), 158 G4ThreeVector(0), Mom); 159 160 aTrack->SetPosition(theStrings->operator[](astring)->GetKineticTrack()->GetPosition()); 161 162 #ifdef debug_G4ExcitedStringDecay 163 G4cout<<" A particle stored in the track is "<<aTrack->GetDefinition()->GetParticleName()<<G4endl; 164 #endif 165 166 generatedKineticTracks = new G4KineticTrackVector; 167 generatedKineticTracks->push_back(aTrack); 168 } 169 170 if (generatedKineticTracks == nullptr || generatedKineticTracks->size() == 0) 171 { 172 // G4cerr << "G4VPartonStringModel:No KineticTracks produced" << G4endl; 173 delete generatedKineticTracks; 174 generatedKineticTracks = nullptr; 175 continue; 176 } 177 178 G4LorentzVector KTsum1(0.,0.,0.,0.); 179 for ( unsigned int aTrack=0; aTrack<generatedKineticTracks->size();aTrack++) 180 { 181 #ifdef debug_G4ExcitedStringDecay 182 G4cout<<"Prod part No. "<<aTrack+1<<" " 183 <<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" " 184 <<(*generatedKineticTracks)[aTrack]->Get4Momentum() 185 <<(*generatedKineticTracks)[aTrack]->Get4Momentum().mag()<<G4endl; 186 #endif 187 // --------------- Sampling mass of unstable hadronic resonances ---------------- 188 TrackDefinition = (*generatedKineticTracks)[aTrack]->GetDefinition(); 189 190 if (TrackDefinition->IsShortLived()) 191 { 192 G4double NewTrackMass = 193 BrW.SampleMass( TrackDefinition->GetPDGMass(), TrackDefinition->GetPDGWidth(), 194 BrW.GetMinimumMass( TrackDefinition ) + 10.0*MeV, 195 TrackDefinition->GetPDGMass() + 5.0*TrackDefinition->GetPDGWidth() ); 196 G4LorentzVector Tmp=G4LorentzVector((*generatedKineticTracks)[aTrack]->Get4Momentum()); 197 Tmp.setE(std::sqrt(sqr(NewTrackMass) + Tmp.vect().mag2())); 198 199 (*generatedKineticTracks)[aTrack]->Set4Momentum(Tmp); 200 201 #ifdef debug_G4ExcitedStringDecay 202 G4cout<<"Resonance *** "<<aTrack+1<<" " 203 <<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" " 204 <<(*generatedKineticTracks)[aTrack]->Get4Momentum() 205 <<(*generatedKineticTracks)[aTrack]->Get4Momentum().mag()<<G4endl; 206 #endif 207 } 208 //------------------------------------------------------------------------------- 209 210 theResult->push_back(generatedKineticTracks->operator[](aTrack)); 211 KTsum1+= (*generatedKineticTracks)[aTrack]->Get4Momentum(); 212 } 213 KTsecondaries+=KTsum1; 214 215 #ifdef debug_G4ExcitedStringDecay 216 G4cout << "String secondaries(" <<generatedKineticTracks->size()<< ")"<<G4endl 217 <<"Init string momentum: "<< theStrings->operator[](astring)->Get4Momentum()<<G4endl 218 <<"Final hadrons momentum: "<< KTsum1 << G4endl; 219 #endif 220 221 if ( KTsum1.e() > 0 && 222 std::abs((KTsum1.e()-theStrings->operator[](astring)->Get4Momentum().e()) / KTsum1.e()) > perMillion ) 223 { 224 NeedEnergyCorrector=true; 225 } 226 227 #ifdef debug_G4ExcitedStringDecay 228 G4cout<<"NeedEnergyCorrection yes/no "<<NeedEnergyCorrector<<G4endl; 229 #endif 230 231 // clean up 232 delete generatedKineticTracks; 233 success=true; 234 } 235 236 if ( NeedEnergyCorrector ) success=EnergyAndMomentumCorrector(theResult, KTsum); 237 } while (!success && (attempts < 100)); /* Loop checking, 07.08.2015, A.Ribon */ 238 239 for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++) 240 { 241 Ptmp=(*theResult)[aTrack]->Get4Momentum(); 242 Ptmp.transform( toLab); 243 (*theResult)[aTrack]->Set4Momentum(Ptmp); 244 } 245 246 #ifdef debug_G4ExcitedStringDecay 247 G4cout<<"End of the strings fragmentation (G4ExcitedStringDecay)"<<G4endl; 248 249 G4LorentzVector KTsum1(0.,0.,0.,0.); 250 251 for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++) 252 { 253 G4cout << " corrected tracks .. " << (*theResult)[aTrack]->GetDefinition()->GetParticleName() 254 <<" " << (*theResult)[aTrack]->Get4Momentum() 255 <<" " << (*theResult)[aTrack]->Get4Momentum().mag()<< G4endl; 256 KTsum1+= (*theResult)[aTrack]->Get4Momentum(); 257 } 258 259 G4cout << "Needcorrector/success " << NeedEnergyCorrector << "/" << success 260 << ", Corrected total 4 momentum " << KTsum1 << G4endl; 261 if ( ! success ) G4cout << "failed to correct E/p" << G4endl; 262 263 G4cout<<"End of the Hadronization (G4ExcitedStringDecay)"<<G4endl; 264 #endif 265 266 if (!success) 267 { 268 if (theResult->size() != 0) 269 { 270 std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack()); 271 theResult->clear(); 272 delete theResult; theResult=0; 273 } 274 for ( unsigned int astring=0; astring < theStrings->size(); astring++) 275 // for ( unsigned int astring=0; astring < 1; astring++) // Need more correct. For testing purposes. 276 { 277 if ( theStrings->operator[](astring)->IsExcited() ) 278 { 279 Ptmp=theStrings->operator[](astring)->GetLeftParton()->Get4Momentum(); 280 Ptmp.transform( toLab); 281 theStrings->operator[](astring)->GetLeftParton()->Set4Momentum(Ptmp); 282 283 Ptmp=theStrings->operator[](astring)->GetRightParton()->Get4Momentum(); 284 Ptmp.transform( toLab); 285 theStrings->operator[](astring)->GetRightParton()->Set4Momentum(Ptmp); 286 } 287 else 288 { 289 Ptmp=theStrings->operator[](astring)->GetKineticTrack()->Get4Momentum(); 290 Ptmp.transform( toLab); 291 theStrings->operator[](astring)->GetKineticTrack()->Set4Momentum(Ptmp); 292 } 293 } 294 } 295 return theResult; 296 } 297 298 299 G4bool G4ExcitedStringDecay:: 300 EnergyAndMomentumCorrector(G4KineticTrackVector* Output, G4LorentzVector& TotalCollisionMom) 301 { 302 const int nAttemptScale = 500; 303 const double ErrLimit = 1.E-5; 304 if (Output->empty()) return TRUE; 305 G4LorentzVector SumMom; 306 G4double SumMass = 0; 307 G4double TotalCollisionMass = TotalCollisionMom.m(); 308 309 std::vector<G4double> HadronMass; G4double HadronM(0.); 310 311 #ifdef debug_G4ExcitedStringCorr 312 G4cout<<G4endl<<"EnergyAndMomentumCorrector. Number of particles: "<<Output->size()<<G4endl; 313 #endif 314 // Calculate sum hadron 4-momenta and summing hadron mass 315 unsigned int cHadron; 316 for (cHadron = 0; cHadron < Output->size(); cHadron++) 317 { 318 SumMom += Output->operator[](cHadron)->Get4Momentum(); 319 HadronM=Output->operator[](cHadron)->Get4Momentum().mag(); HadronMass.push_back(HadronM); 320 SumMass += Output->operator[](cHadron)->Get4Momentum().mag(); //GetDefinition()->GetPDGMass(); 321 } 322 323 #ifdef debug_G4ExcitedStringCorr 324 G4cout<<"Sum part mom "<<SumMom<<" "<<SumMom.mag()<<G4endl 325 <<"Sum str mom "<<TotalCollisionMom<<" "<<TotalCollisionMom.mag()<<G4endl; 326 G4cout<<"SumMass TotalCollisionMass "<<SumMass<<" "<<TotalCollisionMass<<G4endl; 327 #endif 328 329 // Cannot correct a single particle 330 if (Output->size() < 2) return FALSE; 331 332 if (SumMass > TotalCollisionMass) return FALSE; 333 SumMass = SumMom.m2(); 334 if (SumMass < 0) return FALSE; 335 SumMass = std::sqrt(SumMass); 336 337 // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s. 338 // G4ThreeVector Beta = -SumMom.boostVector(); 339 G4ThreeVector Beta = -TotalCollisionMom.boostVector(); 340 Output->Boost(Beta); 341 342 // Scale total c.m.s. hadron energy (hadron system mass). 343 // It should be equal interaction mass 344 G4double Scale = 1; 345 G4int cAttempt = 0; 346 G4double Sum = 0; 347 G4bool success = false; 348 for (cAttempt = 0; cAttempt < nAttemptScale; cAttempt++) 349 { 350 Sum = 0; 351 for (cHadron = 0; cHadron < Output->size(); cHadron++) 352 { 353 HadronM = HadronMass.at(cHadron); 354 G4LorentzVector HadronMom = Output->operator[](cHadron)->Get4Momentum(); 355 HadronMom.setVect(Scale*HadronMom.vect()); 356 G4double E = std::sqrt(HadronMom.vect().mag2() + sqr(HadronM)); 357 //sqr(Output->operator[](cHadron)->GetDefinition()->GetPDGMass())); 358 HadronMom.setE(E); 359 Output->operator[](cHadron)->Set4Momentum(HadronMom); 360 Sum += E; 361 } 362 Scale = TotalCollisionMass/Sum; 363 #ifdef debug_G4ExcitedStringCorr 364 G4cout << "Scale-1=" << Scale -1 365 << ", TotalCollisionMass=" << TotalCollisionMass 366 << ", Sum=" << Sum 367 << G4endl; 368 #endif 369 if (std::fabs(Scale - 1) <= ErrLimit) 370 { 371 success = true; 372 break; 373 } 374 } 375 376 #ifdef debug_G4ExcitedStringCorr 377 if (!success) 378 { 379 G4cout << "G4ExcitedStringDecay::EnergyAndMomentumCorrector - Warning"<<G4endl; 380 G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl; 381 G4cout << " Number of secondaries: " << Output->size() << G4endl; 382 G4cout << " Wanted total energy: " << TotalCollisionMom.e() << G4endl; 383 G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl; 384 // throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay failed to correct..."); 385 } 386 #endif 387 388 // Compute c.m.s. interaction velocity and KTV back boost 389 Beta = TotalCollisionMom.boostVector(); 390 Output->Boost(Beta); 391 392 return success; 393 } 394 395