source: trunk/source/processes/hadronic/models/parton_string/hadronization/src/G4ExcitedStringDecay.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 9.5 KB
Line 
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#include "G4ExcitedStringDecay.hh"
28#include "G4KineticTrack.hh"
29
30G4ExcitedStringDecay::G4ExcitedStringDecay()
31{
32 theStringDecay=NULL;
33}
34
35G4ExcitedStringDecay::G4ExcitedStringDecay(G4VLongitudinalStringDecay * aStringDecay)
36: theStringDecay(aStringDecay)
37{}
38
39G4ExcitedStringDecay::G4ExcitedStringDecay(const G4ExcitedStringDecay &) : G4VStringFragmentation()
40{
41}
42
43G4ExcitedStringDecay::~G4ExcitedStringDecay()
44{
45}
46
47const G4ExcitedStringDecay & G4ExcitedStringDecay::operator=(const G4ExcitedStringDecay &)
48{
49 throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::operator= meant to not be accessable");
50 return *this;
51}
52
53int G4ExcitedStringDecay::operator==(const G4ExcitedStringDecay &) const
54{
55 return 0;
56}
57
58int G4ExcitedStringDecay::operator!=(const G4ExcitedStringDecay &) const
59{
60 return 1;
61}
62
63G4KineticTrackVector *G4ExcitedStringDecay::FragmentString
64 (const G4ExcitedString &theString)
65{
66 if ( theStringDecay == NULL )
67
68 theStringDecay=new G4LundStringFragmentation();
69
70 return theStringDecay->FragmentString(theString);
71}
72
73G4KineticTrackVector *G4ExcitedStringDecay::FragmentStrings
74 (const G4ExcitedStringVector * theStrings)
75{
76 G4LorentzVector KTsum(0.,0.,0.,0.);
77 for ( unsigned int astring=0; astring < theStrings->size(); astring++)
78 {
79 KTsum+= theStrings->operator[](astring)->Get4Momentum();
80
81 if( !(KTsum.e()<1) && !(KTsum.e()>-1) )
82 {
83 throw G4HadronicException(__FILE__, __LINE__,
84 "G4ExcitedStringDecay::FragmentStrings received nan string...");
85 }
86 }
87
88 G4KineticTrackVector * theResult = new G4KineticTrackVector;
89 G4int attempts(0);
90 G4bool success=false;
91 do {
92
93 std::for_each(theResult->begin() , theResult->end() , DeleteKineticTrack());
94 theResult->clear();
95
96 attempts++;
97 G4LorentzVector KTsecondaries(0.,0.,0.,0.);
98 G4bool NeedEnergyCorrector=false;
99
100 for ( unsigned int astring=0; astring < theStrings->size(); astring++)
101 {
102//G4cout<<G4endl<<"String No "<<astring+1<<" "<<theStrings->operator[](astring)->Get4Momentum().mag()<<G4endl;
103
104 G4KineticTrackVector * generatedKineticTracks = NULL;
105
106 if ( theStrings->operator[](astring)->IsExcited() )
107 {
108 generatedKineticTracks=FragmentString(*theStrings->operator[](astring));
109 } else {
110 generatedKineticTracks = new G4KineticTrackVector;
111 generatedKineticTracks->push_back(theStrings->operator[](astring)->GetKineticTrack());
112 }
113
114 if (generatedKineticTracks == NULL)
115 {
116 G4cerr << "G4VPartonStringModel:No KineticTracks produced" << G4endl;
117 continue;
118 }
119
120 G4LorentzVector KTsum1(0.,0.,0.,0.);
121 for ( unsigned int aTrack=0; aTrack<generatedKineticTracks->size();aTrack++)
122 {
123//G4cout<<" Prod part "<<(*generatedKineticTracks)[aTrack]->GetDefinition()->GetParticleName()<<" "<<(*generatedKineticTracks)[aTrack]->Get4Momentum()<<G4endl;
124 theResult->push_back(generatedKineticTracks->operator[](aTrack));
125 KTsum1+= (*generatedKineticTracks)[aTrack]->Get4Momentum();
126 }
127 KTsecondaries+=KTsum1;
128
129 if ( KTsum1.e() > 0 && std::abs((KTsum1.e()-theStrings->operator[](astring)->Get4Momentum().e()) / KTsum1.e()) > perMillion )
130 {
131//--debug-- G4cout << "String secondaries(" <<generatedKineticTracks->size()<< ") momentum: "
132//--debug-- << theStrings->operator[](astring)->Get4Momentum() << " " << KTsum1 << G4endl;
133 NeedEnergyCorrector=true;
134 }
135
136// clean up
137 delete generatedKineticTracks;
138 }
139//--DEBUG G4cout << "Strings/secs total 4 momentum " << KTsum << " " <<KTsecondaries << G4endl;
140
141 success=true;
142 if ( NeedEnergyCorrector ) success=EnergyAndMomentumCorrector(theResult, KTsum);
143
144 } while(!success && (attempts < 100));
145
146#ifdef debug_ExcitedStringDecay
147 G4LorentzVector KTsum1=0;
148 for ( unsigned int aTrack=0; aTrack<theResult->size();aTrack++)
149 {
150 G4cout << " corrected tracks .. " << (*theResult)[aTrack]->GetDefinition()->GetParticleName()
151 <<" " << (*theResult)[aTrack]->Get4Momentum() << G4endl;
152 KTsum1+= (*theResult)[aTrack]->Get4Momentum();
153 }
154 G4cout << "Needcorrector/success " << NeedEnergyCorrector << "/" << success << ", Corrected total 4 momentum " << KTsum1 << G4endl;
155 if ( ! success ) G4cout << "failed to correct E/p" << G4endl;
156#endif
157
158 return theResult;
159}
160
161G4bool G4ExcitedStringDecay::EnergyAndMomentumCorrector
162 (G4KineticTrackVector* Output, G4LorentzVector& TotalCollisionMom)
163 {
164 const int nAttemptScale = 500;
165 const double ErrLimit = 1.E-5;
166 if (Output->empty())
167 return TRUE;
168 G4LorentzVector SumMom;
169 G4double SumMass = 0;
170 G4double TotalCollisionMass = TotalCollisionMom.m();
171
172 if( !(TotalCollisionMass<1) && !(TotalCollisionMass>-1) )
173 {
174 std::cout << "TotalCollisionMomentum = "<<TotalCollisionMom<<G4endl;
175 throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay received nan mass...");
176 }
177
178//G4cout<<G4endl<<"EnergyAndMomentumCorrector "<<Output->size()<<G4endl;
179 // Calculate sum hadron 4-momenta and summing hadron mass
180 unsigned int cHadron;
181 for(cHadron = 0; cHadron < Output->size(); cHadron++)
182 {
183 SumMom += Output->operator[](cHadron)->Get4Momentum();
184 if( !(SumMom<1) && !(SumMom>-1) )
185 {
186 throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::EnergyAndMomentumCorrector() received nan momentum...");
187 }
188 SumMass += Output->operator[](cHadron)->GetDefinition()->GetPDGMass();
189 if( !(SumMass<1) && !(SumMass>-1) )
190 {
191 throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay::EnergyAndMomentumCorrector() received nan mass...");
192 }
193 }
194
195//G4cout<<"SumMass TotalCollisionMass "<<SumMass<<" "<<TotalCollisionMass<<G4endl;
196
197 // Cannot correct a single particle
198 if (Output->size() < 2) return FALSE;
199
200 if (SumMass > TotalCollisionMass) return FALSE;
201 SumMass = SumMom.m2();
202 if (SumMass < 0) return FALSE;
203 SumMass = std::sqrt(SumMass);
204
205 // Compute c.m.s. hadron velocity and boost KTV to hadron c.m.s.
206 G4ThreeVector Beta = -SumMom.boostVector();
207 Output->Boost(Beta);
208
209 // Scale total c.m.s. hadron energy (hadron system mass).
210 // It should be equal interaction mass
211 G4double Scale = 1;
212 G4int cAttempt = 0;
213 G4double Sum = 0;
214 G4bool success = false;
215 for(cAttempt = 0; cAttempt < nAttemptScale; cAttempt++)
216 {
217 Sum = 0;
218 for(cHadron = 0; cHadron < Output->size(); cHadron++)
219 {
220 G4LorentzVector HadronMom = Output->operator[](cHadron)->Get4Momentum();
221 HadronMom.setVect(Scale*HadronMom.vect());
222 G4double E = std::sqrt(HadronMom.vect().mag2() + sqr(Output->operator[](cHadron)->GetDefinition()->GetPDGMass()));
223 HadronMom.setE(E);
224 Output->operator[](cHadron)->Set4Momentum(HadronMom);
225 Sum += E;
226 }
227 Scale = TotalCollisionMass/Sum;
228#ifdef debug_G4ExcitedStringDecay
229 G4cout << "Scale-1=" << Scale -1
230 << ", TotalCollisionMass=" << TotalCollisionMass
231 << ", Sum=" << Sum
232 << G4endl;
233#endif
234 if (std::fabs(Scale - 1) <= ErrLimit)
235 {
236 success = true;
237 break;
238 }
239 }
240#ifdef debug_G4ExcitedStringDecay
241 if(!success)
242 {
243 G4cout << "G4ExcitedStringDecay::EnergyAndMomentumCorrector - Warning"<<G4endl;
244 G4cout << " Scale not unity at end of iteration loop: "<<TotalCollisionMass<<" "<<Sum<<" "<<Scale<<G4endl;
245 G4cout << " Number of secondaries: " << Output->size() << G4endl;
246 G4cout << " Wanted total energy: " << TotalCollisionMom.e() << G4endl;
247 G4cout << " Increase number of attempts or increase ERRLIMIT"<<G4endl;
248// throw G4HadronicException(__FILE__, __LINE__, "G4ExcitedStringDecay failed to correct...");
249 }
250#endif
251 // Compute c.m.s. interaction velocity and KTV back boost
252 Beta = TotalCollisionMom.boostVector();
253 Output->Boost(Beta);
254
255 return success;
256 }
Note: See TracBrowser for help on using the repository browser.