Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QGSParticipants.cc
Go to the documentation of this file.
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#include <utility>
27
28#include "G4QGSParticipants.hh"
29#include "globals.hh"
30#include "G4SystemOfUnits.hh"
31#include "G4LorentzVector.hh"
32#include "G4V3DNucleus.hh"
33#include "G4ParticleTable.hh"
34#include "G4IonTable.hh"
36
37#include "G4Exp.hh"
38#include "G4Log.hh"
39#include "G4Pow.hh"
40
41//#define debugQGSParticipants
42//#define debugPutOnMassShell
43
44// Class G4QGSParticipants
45
46// Promoting model parameters from local variables class properties
48
50 ModelMode(SOFT),
51 nCutMax(7),ThresholdParameter(0.000*GeV),
52 QGSMThreshold(3*GeV),theNucleonRadius(1.5*fermi),alpha(-0.5),beta(2.5)
53{
54 // Parameters setting
55 SetCofNuclearDestruction(1.);
56 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
57 SetDofNuclearDestruction( 0.3 );
58 SetPt2ofNuclearDestruction( 0.075*GeV*GeV );
59 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
60 SetExcitationEnergyPerWoundedNucleon( 40.0*MeV );
61
62 sigmaPt=0.25*sqr(GeV);
63}
64
66: G4VParticipants(),ModelMode(right.ModelMode), nCutMax(right.nCutMax),
67 ThresholdParameter(right.ThresholdParameter), QGSMThreshold(right.QGSMThreshold),
68 theNucleonRadius(right.theNucleonRadius)
69{}
70
72
74{
75 theProjectile = thePrimary;
76
77 Regge = new G4Reggeons(theProjectile.GetDefinition());
78
80
81 NumberOfInvolvedNucleonsOfProjectile= 0;
82 G4LorentzVector tmp( 0.0, 0.0, 0.0, 0.0 );
83 ProjectileResidualMassNumber = 0;
84 ProjectileResidualCharge = 0;
85 ProjectileResidualExcitationEnergy = 0.0;
86 ProjectileResidual4Momentum = tmp;
87
88 NumberOfInvolvedNucleonsOfTarget= 0;
89 TargetResidualMassNumber = theNucleus->GetMassNumber();
90 TargetResidualCharge = theNucleus->GetCharge();
91 TargetResidualExcitationEnergy = 0.0;
92
94 G4Nucleon* NuclearNucleon;
95 while ( ( NuclearNucleon = theNucleus->GetNextNucleon() ) ) {
96 tmp+=NuclearNucleon->Get4Momentum();
97 }
98 TargetResidual4Momentum = tmp;
99
100 if ( std::abs( theProjectile.GetDefinition()->GetBaryonNumber() ) <= 1 ) {
101 // Projectile is a hadron : meson or baryon
102 ProjectileResidualMassNumber = std::abs( theProjectile.GetDefinition()->GetBaryonNumber() );
103 ProjectileResidualCharge = G4int( theProjectile.GetDefinition()->GetPDGCharge() );
104 ProjectileResidualExcitationEnergy = 0.0;
105 ProjectileResidual4Momentum.setVect( theProjectile.GetMomentum() );
106 ProjectileResidual4Momentum.setE( theProjectile.GetTotalEnergy() );
107 }
108
109 #ifdef debugQGSParticipants
110 G4cout <<G4endl<< "G4QGSParticipants::BuildInteractions---------" << G4endl
111 << "thePrimary " << thePrimary.GetDefinition()->GetParticleName()<<" "
112 <<ProjectileResidual4Momentum<<G4endl;
113 G4cout << "Target A and Z " << theNucleus->GetMassNumber() <<" "<<theNucleus->GetCharge()<<" "
114 << TargetResidual4Momentum<<G4endl;
115 #endif
116
117 G4bool Success( true );
118
119 const G4int maxNumberOfLoops = 1000;
120 G4int loopCounter = 0;
121 do
122 {
123 const G4int maxNumberOfInternalLoops = 1000;
124 G4int internalLoopCounter = 0;
125 do
126 {
127 if(std::abs(theProjectile.GetDefinition()->GetPDGEncoding()) < 100.0)
128 {
129 SelectInteractions(theProjectile); // for lepton projectile
130 }
131 else
132 {
133 GetList( theProjectile ); // Get list of participating nucleons for hadron projectile
134 }
135
136 if ( theInteractions.size() == 0 ) return;
137
138 StoreInvolvedNucleon(); // Store participating nucleon
139
140 ReggeonCascade(); // Make reggeon cascading. Involve nucleons in the cascading.
141
142 Success = PutOnMassShell(); // Ascribe momenta to the involved and participating nucleons
143
144 if(!Success) PrepareInitialState( thePrimary );
145
146 } while( (!Success) && ++internalLoopCounter < maxNumberOfInternalLoops );
147
148 if ( internalLoopCounter >= maxNumberOfInternalLoops ) {
149 Success = false;
150 }
151
152 if ( Success ) {
153 #ifdef debugQGSParticipants
154 G4cout<<G4endl<<"PerformDiffractiveCollisions(); if they happend." <<G4endl;
155 #endif
156
158
159 #ifdef debugQGSParticipants
160 G4cout<<G4endl<<"SplitHadrons();" <<G4endl;
161 #endif
162
163 SplitHadrons();
164
166 #ifdef debugQGSParticipants
167 G4cout<<"Perform non-Diffractive Collisions if they happend. Determine Parton Momenta and so on." <<G4endl;
168 #endif
169 Success = DeterminePartonMomenta();
170 }
171
172 if(!Success) PrepareInitialState( thePrimary );
173 }
174 } while( (!Success) && ++loopCounter < maxNumberOfLoops );
175
176 if ( loopCounter >= maxNumberOfLoops ) {
177 Success = false;
178 #ifdef debugQGSParticipants
179 G4cout<<"NOT Successful ======" <<G4endl;
180 #endif
181 }
182
183 if ( Success ) {
184 CreateStrings(); // To create strings
185
186 GetResiduals(); // To calculate residual nucleus properties
187
188 #ifdef debugQGSParticipants
189 G4cout<<"All O.K. ======" <<G4endl;
190 #endif
191 }
192
193 // clean-up, if necessary
194 #ifdef debugQGSParticipants
195 G4cout<<"Clearing "<<G4endl;
196 G4cout<<"theInteractions.size() "<<theInteractions.size()<<G4endl;
197 #endif
198
199 if( Regge ) delete Regge;
200
201 std::for_each( theInteractions.begin(), theInteractions.end(), DeleteInteractionContent() );
202 theInteractions.clear();
203
204 // Erasing of target involved nucleons.
205 #ifdef debugQGSParticipants
206 G4cout<<"Erasing of involved target nucleons "<<NumberOfInvolvedNucleonsOfTarget<<G4endl;
207 #endif
208
209 if ( NumberOfInvolvedNucleonsOfTarget != 0 ) {
210 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; i++ ) {
211 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfTarget[i]->GetSplitableHadron();
212 if ( (aNucleon != 0 ) && (aNucleon->GetStatus() >= 1) ) delete aNucleon;
213 }
214 }
215
216 // Erasing of projectile involved nucleons.
217 if ( NumberOfInvolvedNucleonsOfProjectile != 0 ) {
218 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfProjectile; i++ ) {
219 G4VSplitableHadron* aNucleon = TheInvolvedNucleonsOfProjectile[i]->GetSplitableHadron();
220 if ( aNucleon ) delete aNucleon;
221 }
222 }
223
224 #ifdef debugQGSParticipants
225 G4cout<<"Delition of target nucleons from soft interactions "<<theTargets.size()
226 <<G4endl<<G4endl;
227 #endif
228 std::for_each(theTargets.begin(), theTargets.end(), DeleteSplitableHadron());
229 theTargets.clear();
230
234 }
235}
236
237//===========================================================
238void G4QGSParticipants::PrepareInitialState( const G4ReactionProduct& thePrimary )
239{
240 // Clearing of the arrays
241 // Erasing of the projectile
242 G4InteractionContent* anIniteraction = theInteractions[0];
243 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
244 if( pProjectile ) delete pProjectile;
245
246 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
247 theInteractions.clear();
248
249 // Erasing of the envolved nucleons and target nucleons from diffraction dissociations
251 G4Nucleon* aNucleon;
252 while ( ( aNucleon = theNucleus->GetNextNucleon() ) )
253 {
254 if ( aNucleon->AreYouHit() ) {
255 G4VSplitableHadron* splaNucleon = aNucleon->GetSplitableHadron();
256 if ( (splaNucleon != 0) && (splaNucleon->GetStatus() >=1) ) delete splaNucleon;
257 aNucleon->Hit(nullptr);
258 NumberOfInvolvedNucleonsOfTarget--;
259 }
260 }
261
262 // Erasing of nuclear nucleons participated in soft interactions
263 std::for_each(theTargets.begin(), theTargets.end(), DeleteSplitableHadron());
264 theTargets.clear();
265
266 // Preparation to a new attempt
267 theProjectile = thePrimary;
268
271 DoLorentzBoost(-theCurrentVelocity); // Lorentz boost of the target nucleus
272
273 if (theNucleus->GetMassNumber() == 1)
274 {
275 G4ThreeVector aPos = G4ThreeVector(0.,0.,0.);
278 tNucleon->SetPosition(aPos);
279 }
280
281 G4LorentzVector Tmp( 0.0, 0.0, 0.0, 0.0 );
282 NumberOfInvolvedNucleonsOfTarget= 0;
283 TargetResidualMassNumber = theNucleus->GetMassNumber();
284 TargetResidualCharge = theNucleus->GetCharge();
285 TargetResidualExcitationEnergy = 0.0;
286
287 G4Nucleon* NuclearNucleon;
288 while ( ( NuclearNucleon = theNucleus->GetNextNucleon() ) )
289 {Tmp+=NuclearNucleon->Get4Momentum();}
290
291 TargetResidual4Momentum = Tmp;
292}
293
294//===========================================================
295void G4QGSParticipants::GetList( const G4ReactionProduct& thePrimary ) {
296 #ifdef debugQGSParticipants
297 G4cout<<G4endl<<"G4QGSParticipants::GetList +++++++++++++"<<G4endl;
298 #endif
299
300 // Direction: True - Proj, False - Target
303
304 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
305 G4LorentzVector aNucleonMomentum(0.,0.,0., 938.0*MeV);
306
307 G4double SS=(aPrimaryMomentum + aNucleonMomentum).mag2();
308
309 Regge->SetS(SS);
310
311 //--------------------------------------
313 G4Nucleon * tNucleon = theNucleus->GetNextNucleon();
314
315 if ( ! tNucleon ) {
316 #ifdef debugQGSParticipants
317 G4cout << "QGSM - BAD situation: pNucleon is NULL ! Leaving immediately!" << G4endl;
318 #endif
319 return;
320 }
321
322 G4double theNucleusOuterR = theNucleus->GetOuterRadius();
323
324 if (theNucleus->GetMassNumber() == 1)
325 {
326 G4ThreeVector aPos = G4ThreeVector(0.,0.,0.);
327 tNucleon->SetPosition(aPos);
328 theNucleusOuterR = 0.;
329 }
330
331 // Determination of participating nucleons of nucleus ------------------------------------
332
333 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
334 theInteractions.clear();
335
336 G4int totalCuts = 0;
337 G4int MaxPower=thePrimary.GetMomentum().mag()/(3.3*GeV); if(MaxPower < 1) MaxPower=1;
338
339 const G4int maxNumberOfLoops = 1000;
340
341 G4int NumberOfTries = 0;
342 G4double Scale = 1.0;
343
344 G4int loopCounter = -1;
345 while( (theInteractions.size() == 0) && ++loopCounter < maxNumberOfLoops )
346 {
347 InteractionMode = ALL; // Mode = ALL, WITHOUT_R, NON_DIFF
348
349 // choose random impact parameter of a collision
350 std::pair<G4double, G4double> theImpactParameter;
351
352 NumberOfTries++;
353 if( NumberOfTries == 100*(NumberOfTries/100) ) Scale /=2.0;
354
355 theImpactParameter = theNucleus->ChooseImpactXandY(theNucleusOuterR/Scale + theNucleonRadius);
356 G4double impactX = theImpactParameter.first;
357 G4double impactY = theImpactParameter.second;
358
359 #ifdef debugQGSParticipants
360 G4cout<<"InteractionMode "<<InteractionMode<<G4endl;
361 G4cout<<"Impact parameter (fm ) "<<std::sqrt(sqr(impactX)+sqr(impactY))/fermi<<" "<<G4endl;
362 #endif
363
364 // loop over nucleons to find collisions
366 G4int nucleonCount = -1;
368
369 G4double Power=MaxPower;
370
371 while( (tNucleon = theNucleus->GetNextNucleon()) )
372 {
373 if(Power <= 0.) break;
374 nucleonCount++;
375
376 G4LorentzVector nucleonMomentum=tNucleon->Get4Momentum();
377
378 G4double Distance2 = sqr(impactX - tNucleon->GetPosition().x()) +
379 sqr(impactY - tNucleon->GetPosition().y());
380
381 G4double Pint(0.); // A probability of interaction at given impact parameter
382 G4double Pprd(0.), Ptrd(0.), Pdd(0.); // Probabilities of Proj. diffr., Target diffr., Double diffr.
383 G4double Pnd (0.), Pnvr(0.); // Probabilities of non-diffr. and quark exchange
384 G4int NcutPomerons(0); // Number of cutted pomerons
385
386 Regge->GetProbabilities(std::sqrt(Distance2), InteractionMode,
387 Pint, Pprd, Ptrd, Pdd, Pnd, Pnvr);
388 #ifdef debugQGSParticipants
389 G4cout<<"Nucleon & its impact parameter: "<<nucleonCount<<" "<<std::sqrt(Distance2)/fermi<<" (fm)"<<G4endl;
390 G4cout<<"Probability of interaction: "<<Pint<<G4endl;
391 G4cout<<"Probability of PrD, TrD, DD: "<<Pprd<<" "<<Ptrd<<" "<<Pdd<<G4endl;
392 G4cout<<"Probability of NonDiff, QuarkExc.: "<<Pnd<<" "<<Pnvr<<" in inel. inter."<<G4endl;
393 #endif
394
395 if (Pint > G4UniformRand())
396 { // An interaction is happend.
397
398 G4double rndNumber = G4UniformRand();
399 G4int InteractionType(0);
400
401 if((InteractionMode==ALL)||(InteractionMode==WITHOUT_R)) // Mode = ALL, WITHOUT_R, NON_DIFF
402 {
403 if( rndNumber < Pprd ) {InteractionType = PrD; InteractionMode = WITHOUT_R;}
404 else if( rndNumber < Pprd+Ptrd) {InteractionType = TrD; InteractionMode = WITHOUT_R;}
405 else if( rndNumber < Pprd+Ptrd+Pdd) {InteractionType = DD; InteractionMode = WITHOUT_R;}
406 else if( rndNumber < Pprd+Ptrd+Pdd+Pnd ) {InteractionType = NonD; InteractionMode = NON_DIFF;
407 NcutPomerons = Regge->ncPomerons(); }
408 else {InteractionType = Qexc; InteractionMode = ALL; }
409 }
410 else // InteractionMode == NON_DIFF
411 {
412 InteractionMode = NON_DIFF;
413 if( rndNumber < Ptrd ) {InteractionType = TrD; }
414 else if( rndNumber < Ptrd + Pnd) {InteractionType = NonD; NcutPomerons = Regge->ncPomerons();}
415 }
416
417 if( (InteractionType == NonD) && (NcutPomerons == 0)) continue;
418
420 G4QGSMSplitableHadron* aTargetSPB = new G4QGSMSplitableHadron(*tNucleon);
421 tNucleon->Hit(aTargetSPB);
422
423 #ifdef debugQGSParticipants
424 G4cout<<"An interaction is happend."<<G4endl;
425 G4cout<<"Target nucleon - "<<nucleonCount<<" "
426 <<tNucleon->GetDefinition()->GetParticleName()<<G4endl;
427 G4cout<<"Interaction type:"<<InteractionType
428 <<" (0 -PrD, 1 - TrD, 2 - DD, 3 - NonD, 4 - Qexc)"<<G4endl;
429 G4cout<<"New Inter. mode:"<<InteractionMode
430 <<" (0 -ALL, 1 - WITHOUT_R, 2 - NON_DIFF)"<<G4endl;
431 if( InteractionType == NonD )
432 G4cout<<"Number of cutted pomerons: "<<NcutPomerons<<G4endl;
433 #endif
434
435 if((InteractionType == PrD) || (InteractionType == TrD) || (InteractionType == DD) ||
436 (InteractionType == Qexc))
437 { // diffractive-like interaction occurs
438 #ifdef debugQGSParticipants
439 G4cout<<"Diffractive-like interaction occurs"<<G4endl;
440 #endif
441
444
445 aInteraction->SetTarget(aTargetSPB);
446 aInteraction->SetTargetNucleon(tNucleon);
447 aTargetSPB->SetCollisionCount(0);
448 aTargetSPB->SetStatus(1);
449
450 aInteraction->SetNumberOfDiffractiveCollisions(1);
451 aInteraction->SetNumberOfSoftCollisions(0);
452 aInteraction->SetStatus(InteractionType);
453 theInteractions.push_back(aInteraction);
454 }
455 else
456 { // nondiffractive interaction occurs
457 #ifdef debugQGSParticipants
458 G4cout<<"Non-diffractive interaction occurs, max NcutPomerons "<<NcutPomerons<<G4endl;
459 #endif
460
461 G4int nCuts;
462
463 G4int Vncut=0;
464 for(nCuts = 0; nCuts < NcutPomerons; nCuts++)
465 {
466 if( G4UniformRand() < Power/MaxPower ){Vncut++; Power--; if(Power <= 0.) break;}
467 }
468 nCuts=Vncut;
469
470 if( nCuts == 0 ) {delete aTargetSPB; tNucleon->Hit(nullptr); continue;}
471
472 totalCuts += nCuts;
473 #ifdef debugQGSParticipants
474 G4cout<<"Number of cuts in the interaction "<<nCuts<<G4endl;
475 #endif
476
477 aTargetSPB->IncrementCollisionCount(nCuts);
478 aTargetSPB->SetStatus(0);
479 theTargets.push_back(aTargetSPB);
480
483
485 aInteraction->SetTarget(aTargetSPB);
486 aInteraction->SetTargetNucleon(tNucleon);
487 aInteraction->SetNumberOfSoftCollisions(nCuts);
488 aInteraction->SetStatus(InteractionType);
489 theInteractions.push_back(aInteraction);
490 }
491 } // End of if (Pint > G4UniformRand())
492 } // End of while( (tNucleon = theNucleus->GetNextNucleon()) )
493
494 #ifdef debugQGSParticipants
495 G4cout << G4endl<<"Number of wounded nucleons "<<G4QGSParticipants_NPart<<G4endl;
496 #endif
497
498 } // End of while( (theInteractions.size() == 0) && ++loopCounter < maxNumberOfLoops )
499
500 if ( loopCounter >= maxNumberOfLoops ) {
501 #ifdef debugQGSParticipants
502 G4cout <<"BAD situation: forced loop exit!" << G4endl;
503 #endif
504 // Perhaps there is something to set here...
505 // Decrease impact parameter ??
506 // Select collisions with only diffraction ??
507 // Selecy only non-diffractive interactions ??
508 }
509 //------------------------------------------------------------
510 std::vector<G4InteractionContent*>::iterator i;
511
512 if( theInteractions.size() != 0)
513 {
514 if( InteractionMode == ALL ) // It can be if all interactions were quark-exchange.
515 { // Only the first one will be saved, all other will be erased.
516 i = theInteractions.end()-1;
517
518 while ( theInteractions.size() != 1 )
519 {
520 G4InteractionContent* anInteraction = *i;
521 G4Nucleon * pNucleon = anInteraction->GetTargetNucleon(); pNucleon->Hit(nullptr);
522 delete anInteraction->GetTarget();
523 delete *i;
524 i=theInteractions.erase(i);
525 i--;
526 }
527 }
528 else
529 { // All quark exchanges will be erased
530 i = theInteractions.begin();
531 while ( i != theInteractions.end() )
532 {
533 G4InteractionContent* anInteraction = *i;
534
535 if( anInteraction->GetStatus() == Qexc )
536 {
537 G4Nucleon* aTargetNucleon = anInteraction->GetTargetNucleon();
538 aTargetNucleon->Hit(nullptr);
539
540 delete anInteraction->GetTarget();
541 delete *i;
542 i=theInteractions.erase(i);
543 }
544 else
545 {
546 i++;
547 }
548 }
549 }
550
551 #ifdef debugQGSParticipants
552 G4cout <<"Total number of cuts "<< totalCuts <<G4endl;
553 #endif
554 }
555}
556
557//=============================================================
558void G4QGSParticipants::StoreInvolvedNucleon()
559{ //To store nucleons involved in the interaction
560
561 NumberOfInvolvedNucleonsOfTarget = 0;
562
564
565 G4Nucleon* aNucleon;
566 while ( ( aNucleon = theNucleus->GetNextNucleon() ) ) {
567 if ( aNucleon->AreYouHit() ) {
568 TheInvolvedNucleonsOfTarget[NumberOfInvolvedNucleonsOfTarget] = aNucleon;
569 NumberOfInvolvedNucleonsOfTarget++;
570 }
571 }
572
573 #ifdef debugQGSParticipants
574 G4cout << G4endl<<"G4QGSParticipants::StoreInvolvedNucleon() if they were "<<G4endl
575 <<"Stored # of wounded nucleons of target "
576 << NumberOfInvolvedNucleonsOfTarget <<G4endl;
577 #endif
578 return;
579}
580
581//=============================================================
582
583void G4QGSParticipants::ReggeonCascade()
584{ // Implementation of the reggeon theory inspired model of nuclear destruction
585 #ifdef debugQGSParticipants
586 G4cout << G4endl<<"Reggeon cascading ........."<<G4endl;
587 G4cout<<"C of nucl. desctruction "<<GetCofNuclearDestruction()
588 <<" R2 "<<GetR2ofNuclearDestruction()/fermi/fermi<<" fermi^2"<<G4endl;
589 #endif
590
591 G4int InitNINt = NumberOfInvolvedNucleonsOfTarget;
592
593 // Reggeon cascading in target nucleus
594 for ( G4int InvTN = 0; InvTN < InitNINt; InvTN++ ) {
595 G4Nucleon* aTargetNucleon = TheInvolvedNucleonsOfTarget[ InvTN ];
596
597 G4double CreationTime = aTargetNucleon->GetSplitableHadron()->GetTimeOfCreation();
598
599 G4double XofWoundedNucleon = aTargetNucleon->GetPosition().x();
600 G4double YofWoundedNucleon = aTargetNucleon->GetPosition().y();
601
602 G4V3DNucleus* theTargetNucleus = theNucleus;
603 theTargetNucleus->StartLoop();
604
605 G4int TrgNuc=0;
606 G4Nucleon* Neighbour(0);
607 while ( ( Neighbour = theTargetNucleus->GetNextNucleon() ) ) {
608 TrgNuc++;
609 if ( ! Neighbour->AreYouHit() ) {
610 G4double impact2 = sqr( XofWoundedNucleon - Neighbour->GetPosition().x() ) +
611 sqr( YofWoundedNucleon - Neighbour->GetPosition().y() );
612
613 if ( G4UniformRand() < GetCofNuclearDestruction() *
614 G4Exp( -impact2 / GetR2ofNuclearDestruction() )
615 ) {
616 // The neighbour nucleon is involved in the reggeon cascade
617 #ifdef debugQGSParticipants
618 G4cout<<"Target nucleon involved in reggeon cascading No "<<TrgNuc<<" "
619 <<Neighbour->GetDefinition()->GetParticleName()<<G4endl;
620 #endif
621 TheInvolvedNucleonsOfTarget[ NumberOfInvolvedNucleonsOfTarget ] = Neighbour;
622 NumberOfInvolvedNucleonsOfTarget++;
623
624 G4QGSMSplitableHadron* targetSplitable = new G4QGSMSplitableHadron( *Neighbour );
625
626 Neighbour->Hit( targetSplitable );
627 targetSplitable->SetTimeOfCreation( CreationTime );
628 targetSplitable->SetStatus( 2 );
629 targetSplitable->SetCollisionCount(0);
630
632 anInteraction->SetTarget(targetSplitable);
633 anInteraction->SetTargetNucleon(Neighbour);
634
635 anInteraction->SetNumberOfDiffractiveCollisions(1);
636 anInteraction->SetNumberOfSoftCollisions(0);
637 anInteraction->SetStatus(3);
638 theInteractions.push_back(anInteraction);
639 }
640 }
641 }
642 }
643
644 #ifdef debugQGSParticipants
645 G4cout <<"Number of new involved nucleons "<<NumberOfInvolvedNucleonsOfTarget - InitNINt<<G4endl;
646 #endif
647 return;
648}
649
650//============================================================================
651
652G4bool G4QGSParticipants::PutOnMassShell() {
653
654 G4bool isProjectileNucleus = false;
655 if ( GetProjectileNucleus() ) {
656 isProjectileNucleus = true;
657 }
658
659 #ifdef debugPutOnMassShell
660 G4cout <<G4endl<< "PutOnMassShell start ..............." << G4endl;
661 if ( isProjectileNucleus ) {G4cout << "PutOnMassShell for Nucleus_Nucleus " << G4endl;}
662 #endif
663
664 G4LorentzVector Pprojectile( theProjectile.GetMomentum(), theProjectile.GetTotalEnergy() );
665 if ( Pprojectile.z() < 0.0 ) {
666 return false;
667 }
668
669 G4bool isOk = true;
670
671 G4LorentzVector Ptarget( 0.0, 0.0, 0.0, 0.0 );
672 G4LorentzVector PtargetResidual( 0.0, 0.0, 0.0, 0.0 );
673 G4double SumMasses = 0.0;
674 G4V3DNucleus* theTargetNucleus = GetTargetNucleus();
675 G4double TargetResidualMass = 0.0;
676
677 #ifdef debugPutOnMassShell
678 G4cout << "Target : ";
679 #endif
680
681 isOk = ComputeNucleusProperties( theTargetNucleus, Ptarget, PtargetResidual, SumMasses,
682 TargetResidualExcitationEnergy, TargetResidualMass,
683 TargetResidualMassNumber, TargetResidualCharge );
684
685 if ( ! isOk ) return false;
686
687 G4double Mprojectile = 0.0;
688 G4double M2projectile = 0.0;
689 G4LorentzVector Pproj( 0.0, 0.0, 0.0, 0.0 );
690 G4LorentzVector PprojResidual( 0.0, 0.0, 0.0, 0.0 );
691 G4V3DNucleus* thePrNucleus = GetProjectileNucleus();
692 G4double PrResidualMass = 0.0;
693
694 if ( ! isProjectileNucleus ) { // hadron-nucleus collision
695 Mprojectile = Pprojectile.mag();
696 M2projectile = Pprojectile.mag2();
697 SumMasses += Mprojectile + 20.0*MeV; // Maybe DM must be larger?
698 } else { // nucleus-nucleus or antinucleus-nucleus collision
699
700 #ifdef debugPutOnMassShell
701 G4cout << "Projectile : ";
702 #endif
703
704 isOk = ComputeNucleusProperties( thePrNucleus, Pproj, PprojResidual, SumMasses,
705 ProjectileResidualExcitationEnergy, PrResidualMass,
706 ProjectileResidualMassNumber, ProjectileResidualCharge );
707 if ( ! isOk ) return false;
708 }
709
710 G4LorentzVector Psum = Pprojectile + Ptarget;
711 G4double SqrtS = Psum.mag();
712 G4double S = Psum.mag2();
713
714 #ifdef debugPutOnMassShell
715 G4cout << "Pproj "<<Pprojectile<<G4endl;
716 G4cout << "Ptarg "<<Ptarget<<G4endl;
717 G4cout << "Psum " << Psum/GeV << " GeV" << G4endl << "SqrtS " << SqrtS/GeV << " GeV" << G4endl
718 << "SumMasses, PrResidualMass and TargetResidualMass " << SumMasses/GeV << " "
719 << PrResidualMass/GeV << " " << TargetResidualMass/GeV << " GeV" << G4endl;
720 G4cout << "Ptar res. "<<PtargetResidual<<G4endl;
721 #endif
722
723 if ( SqrtS < SumMasses ) {
724 return false; // It is impossible to simulate after putting nuclear nucleons on mass-shell.
725 }
726
727 // Try to consider also the excitation energy of the residual nucleus, if this is
728 // possible, with the available energy; otherwise, set the excitation energy to zero.
729
730 G4double savedSumMasses = SumMasses;
731 if ( isProjectileNucleus ) {
732 SumMasses -= std::sqrt( sqr( PrResidualMass ) + PprojResidual.perp2() );
733 SumMasses += std::sqrt( sqr( PrResidualMass + ProjectileResidualExcitationEnergy )
734 + PprojResidual.perp2() );
735 }
736 SumMasses -= std::sqrt( sqr( TargetResidualMass ) + PtargetResidual.perp2() );
737 SumMasses += std::sqrt( sqr( TargetResidualMass + TargetResidualExcitationEnergy )
738 + PtargetResidual.perp2() );
739
740 if ( SqrtS < SumMasses ) {
741 SumMasses = savedSumMasses;
742 if ( isProjectileNucleus ) {
743 ProjectileResidualExcitationEnergy = 0.0;
744 }
745 TargetResidualExcitationEnergy = 0.0;
746 }
747
748 TargetResidualMass += TargetResidualExcitationEnergy;
749 if ( isProjectileNucleus ) {
750 PrResidualMass += ProjectileResidualExcitationEnergy;
751 }
752
753 #ifdef debugPutOnMassShell
754 if ( isProjectileNucleus ) {
755 G4cout << "PrResidualMass ProjResidualExcitationEnergy " << PrResidualMass/GeV << " "
756 << ProjectileResidualExcitationEnergy << " MeV" << G4endl;
757 }
758 G4cout << "TargetResidualMass TargetResidualExcitationEnergy " << TargetResidualMass/GeV << " GeV "
759 << TargetResidualExcitationEnergy << " MeV" << G4endl
760 << "Sum masses " << SumMasses/GeV << G4endl;
761 #endif
762
763 // Sampling of nucleons what can transfer to delta-isobars
764 if ( isProjectileNucleus && thePrNucleus->GetMassNumber() != 1 ) {
765 isOk = GenerateDeltaIsobar( SqrtS, NumberOfInvolvedNucleonsOfProjectile,
766 TheInvolvedNucleonsOfProjectile, SumMasses );
767 }
768 if ( theTargetNucleus->GetMassNumber() != 1 ) {
769 isOk = isOk &&
770 GenerateDeltaIsobar( SqrtS, NumberOfInvolvedNucleonsOfTarget,
771 TheInvolvedNucleonsOfTarget, SumMasses );
772 }
773 if ( ! isOk ) return false;
774
775 // Now we know that it is kinematically possible to produce a final state made
776 // of the involved nucleons (or corresponding delta-isobars) and a residual nucleus.
777 // We have to sample the kinematical variables which will allow to define the 4-momenta
778 // of the final state. The sampled kinematical variables refer to the center-of-mass frame.
779 // Notice that the sampling of the transverse momentum corresponds to take into account
780 // Fermi motion.
781
782 // If target is nucleon - return ?
783
784 G4LorentzRotation toCms( -1*Psum.boostVector() );
785 G4LorentzVector Ptmp = toCms*Pprojectile;
786 if ( Ptmp.pz() <= 0.0 ) { // "String" moving backwards in c.m.s., abort collision!
787 return false;
788 }
789
790 G4LorentzRotation toLab( toCms.inverse() );
791
792 G4double YprojectileNucleus = 0.0;
793 if ( isProjectileNucleus ) {
794 Ptmp = toCms*Pproj;
795 YprojectileNucleus = Ptmp.rapidity();
796 }
797 Ptmp = toCms*Ptarget;
798 G4double YtargetNucleus = Ptmp.rapidity();
799
800 // Ascribing of the involved nucleons Pt and Xminus
801 G4double DcorP = 0.0;
802 if ( isProjectileNucleus ) {
803 DcorP = GetDofNuclearDestruction() / thePrNucleus->GetMassNumber();
804 }
805 G4double DcorT = GetDofNuclearDestruction() / theTargetNucleus->GetMassNumber();
806 G4double AveragePt2 = GetPt2ofNuclearDestruction();
807 G4double maxPtSquare = GetMaxPt2ofNuclearDestruction();
808
809 #ifdef debugPutOnMassShell
810 if ( isProjectileNucleus ) {
811 G4cout << "Y projectileNucleus " << YprojectileNucleus << G4endl;
812 }
813 G4cout << "Y targetNucleus " << YtargetNucleus << G4endl
814 << "Dcor " << GetDofNuclearDestruction()
815 << " DcorP DcorT " << DcorP << " " << DcorT << " AveragePt2 " << AveragePt2 << G4endl;
816 #endif
817
818 G4double M2proj = M2projectile; // Initialization needed only for hadron-nucleus collisions
819 G4double WplusProjectile = 0.0;
820 G4double M2target = 0.0;
821 G4double WminusTarget = 0.0;
822 G4int NumberOfTries = 0;
823 G4double ScaleFactor = 1.0;
824 G4bool OuterSuccess = true;
825
826 const G4int maxNumberOfLoops = 1000;
827 G4int loopCounter = 0;
828 do {
829 G4double sqrtM2proj = 0.0, sqrtM2target = 0.0;
830 OuterSuccess = true;
831 const G4int maxNumberOfTries = 1000;
832 do {
833 NumberOfTries++;
834 if ( NumberOfTries == 100*(NumberOfTries/100) ) {
835 // After many tries, it is convenient to reduce the values of DcorP, DcorT and
836 // AveragePt2, so that the sampled momenta (respectively, pz, and pt) of the
837 // involved nucleons (or corresponding delta-isomers) are smaller, and therefore
838 // it is more likely to satisfy the momentum conservation.
839 ScaleFactor /= 2.0;
840 DcorP *= ScaleFactor;
841 DcorT *= ScaleFactor;
842 AveragePt2 *= ScaleFactor;
843 }
844 if ( isProjectileNucleus ) {
845 // Sampling of kinematical properties of projectile nucleons
846 isOk = SamplingNucleonKinematics( AveragePt2, maxPtSquare, DcorP,
847 thePrNucleus, PprojResidual,
848 PrResidualMass, ProjectileResidualMassNumber,
849 NumberOfInvolvedNucleonsOfProjectile,
850 TheInvolvedNucleonsOfProjectile, M2proj );
851 }
852 // Sampling of kinematical properties of target nucleons
853 isOk = isOk &&
854 SamplingNucleonKinematics( AveragePt2, maxPtSquare, DcorT,
855 theTargetNucleus, PtargetResidual,
856 TargetResidualMass, TargetResidualMassNumber,
857 NumberOfInvolvedNucleonsOfTarget,
858 TheInvolvedNucleonsOfTarget, M2target );
859
860 if ( M2proj < 0.0 ) {
862 ed << "Projectile " << theProjectile.GetDefinition()->GetParticleName()
863 << " Target (Z,A)=(" << theTargetNucleus->GetCharge() << "," << theTargetNucleus->GetMassNumber()
864 << ") M2proj=" << M2proj << " -> sets it to 0.0 !" << G4endl;
865 G4Exception( "G4QGSParticipants::PutOnMassShell(): negative projectile squared mass!",
866 "HAD_QGSPARTICIPANTS_002", JustWarning, ed );
867 M2proj = 0.0;
868 }
869 sqrtM2proj = std::sqrt( M2proj );
870 if ( M2target < 0.0 ) {
872 ed << "Projectile " << theProjectile.GetDefinition()->GetParticleName()
873 << " Target (Z,A)=(" << theTargetNucleus->GetCharge() << "," << theTargetNucleus->GetMassNumber()
874 << ") M2target=" << M2target << " -> sets it to 0.0 !" << G4endl;
875 G4Exception( "G4QGSParticipants::PutOnMassShell(): negative target squared mass!",
876 "HAD_QGSPARTICIPANTS_003", JustWarning, ed );
877 M2target = 0.0;
878 };
879 sqrtM2target = std::sqrt( M2target );
880
881 #ifdef debugPutOnMassShell
882 G4cout << "SqrtS, Mp+Mt, Mp, Mt " << SqrtS/GeV << " "
883 << ( sqrtM2proj + sqrtM2target )/GeV << " "
884 << sqrtM2proj/GeV << " " << sqrtM2target/GeV << G4endl;
885 #endif
886
887 if ( ! isOk ) return false;
888 } while ( ( SqrtS < ( sqrtM2proj + sqrtM2target ) ) &&
889 ++NumberOfTries < maxNumberOfTries ); /* Loop checking, 07.08.2015, A.Ribon */
890 if ( NumberOfTries >= maxNumberOfTries ) {
891 return false;
892 }
893 if ( isProjectileNucleus ) {
894 isOk = CheckKinematics( S, SqrtS, M2proj, M2target, YprojectileNucleus, true,
895 NumberOfInvolvedNucleonsOfProjectile,
896 TheInvolvedNucleonsOfProjectile,
897 WminusTarget, WplusProjectile, OuterSuccess );
898 }
899 isOk = isOk &&
900 CheckKinematics( S, SqrtS, M2proj, M2target, YtargetNucleus, false,
901 NumberOfInvolvedNucleonsOfTarget, TheInvolvedNucleonsOfTarget,
902 WminusTarget, WplusProjectile, OuterSuccess );
903 if ( ! isOk ) return false;
904 } while ( ( ! OuterSuccess ) &&
905 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
906 if ( loopCounter >= maxNumberOfLoops ) {
907 return false;
908 }
909
910 // Now the sampling is completed, and we can determine the kinematics of the
911 // whole system. This is done first in the center-of-mass frame, and then it is boosted
912 // to the lab frame. The transverse momentum of the residual nucleus is determined as
913 // the recoil of each hadron (nucleon or delta) which is emitted, i.e. in such a way
914 // to conserve (by construction) the transverse momentum.
915
916 if ( ! isProjectileNucleus ) { // hadron-nucleus collision
917
918 G4double Pzprojectile = WplusProjectile/2.0 - M2projectile/2.0/WplusProjectile;
919 G4double Eprojectile = WplusProjectile/2.0 + M2projectile/2.0/WplusProjectile;
920 Pprojectile.setPz( Pzprojectile );
921 Pprojectile.setE( Eprojectile );
922
923 #ifdef debugPutOnMassShell
924 G4cout << "Proj after in CMS " << Pprojectile/GeV <<" GeV"<< G4endl;
925 #endif
926
927 Pprojectile.transform( toLab );
928 theProjectile.SetMomentum( Pprojectile.vect() );
929 theProjectile.SetTotalEnergy( Pprojectile.e() );
930
932
933 #ifdef debugPutOnMassShell
934 G4cout << "Final proj. mom in Lab. " <<theProjectile.GetMomentum()/GeV<<" "
935 <<theProjectile.GetTotalEnergy()/GeV<<" GeV"<<G4endl;
936 #endif
937
938 } else { // nucleus-nucleus or antinucleus-nucleus collision
939
940 isOk = FinalizeKinematics( WplusProjectile, true, toLab, PrResidualMass,
941 ProjectileResidualMassNumber, NumberOfInvolvedNucleonsOfProjectile,
942 TheInvolvedNucleonsOfProjectile, ProjectileResidual4Momentum );
943
944 #ifdef debugPutOnMassShell
945 G4cout << "Projectile Residual4Momentum in CMS " << ProjectileResidual4Momentum/GeV <<" GeV"<< G4endl;
946 #endif
947
948 if ( ! isOk ) return false;
949
950 ProjectileResidual4Momentum.transform( toLab );
951
952 #ifdef debugPutOnMassShell
953 G4cout << "Projectile Residual4Momentum in Lab " << ProjectileResidual4Momentum/GeV <<" GeV"<< G4endl;
954 #endif
955
956 }
957
958 isOk = FinalizeKinematics( WminusTarget, false, toLab, TargetResidualMass,
959 TargetResidualMassNumber, NumberOfInvolvedNucleonsOfTarget,
960 TheInvolvedNucleonsOfTarget, TargetResidual4Momentum );
961
962 #ifdef debugPutOnMassShell
963 G4cout << "Target Residual4Momentum in CMS " << TargetResidual4Momentum/GeV << " GeV "<< G4endl;
964 #endif
965
966 if ( ! isOk ) return false;
967
968 TargetResidual4Momentum.transform( toLab );
969
970 #ifdef debugPutOnMassShell
971 G4cout << "Target Residual4Momentum in Lab " << TargetResidual4Momentum/GeV << " GeV "<< G4endl;
972 #endif
973
974 return true;
975
976}
977
978//============================================================================
979
980G4ThreeVector G4QGSParticipants::GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const {
981 // @@ this method is used in FTFModel as well. Should go somewhere common!
982
983 G4double Pt2( 0.0 ), Pt(0.0);
984 if ( AveragePt2 > 0.0 ) {
985 G4double x = maxPtSquare/AveragePt2;
986 Pt2 = (x < 200) ?
987 -AveragePt2 * G4Log( 1.0 + G4UniformRand() * ( G4Exp( -x ) -1.0 ) )
988 : -AveragePt2 * G4Log( 1.0 - G4UniformRand() );
989 Pt = std::sqrt( Pt2 );
990 }
991 G4double phi = G4UniformRand() * twopi;
992
993 return G4ThreeVector( Pt*std::cos(phi), Pt*std::sin(phi), 0.0 );
994}
995//============================================================================
996
997G4bool G4QGSParticipants::
998ComputeNucleusProperties( G4V3DNucleus* nucleus, // input parameter
999 G4LorentzVector& nucleusMomentum, // input & output parameter
1000 G4LorentzVector& residualMomentum, // input & output parameter
1001 G4double& sumMasses, // input & output parameter
1002 G4double& residualExcitationEnergy, // input & output parameter
1003 G4double& residualMass, // input & output parameter
1004 G4int& residualMassNumber, // input & output parameter
1005 G4int& residualCharge ) { // input & output parameter
1006
1007 // This method, which is called only by PutOnMassShell, computes some nucleus properties for:
1008 // - either the target nucleus (which is never an antinucleus): this for any kind
1009 // of hadronic interaction (hadron-nucleus, nucleus-nucleus, antinucleus-nucleus);
1010 // - or the projectile nucleus or antinucleus: this only in the case of nucleus-nucleus
1011 // or antinucleus-nucleus interaction.
1012 // This method assumes that the all the parameters have been initialized by the caller;
1013 // the action of this method consists in modifying all these parameters, except the
1014 // first one. The return value is "false" only in the case the pointer to the nucleus
1015 // is null.
1016
1017 if ( ! nucleus ) return false;
1018
1019 G4double ExcitationEPerWoundedNucleon = GetExcitationEnergyPerWoundedNucleon();
1020
1021 // Loop over the nucleons of the nucleus.
1022 // The nucleons that have been involved in the interaction (either from Glauber or
1023 // Reggeon Cascading) will be candidate to be emitted.
1024 // All the remaining nucleons will be the nucleons of the candidate residual nucleus.
1025 // The variable sumMasses is the amount of energy corresponding to:
1026 // 1. transverse mass of each involved nucleon
1027 // 2. 20.0*MeV separation energy for each involved nucleon
1028 // 3. transverse mass of the residual nucleus
1029 // In this first evaluation of sumMasses, the excitation energy of the residual nucleus
1030 // (residualExcitationEnergy, estimated by adding a constant value to each involved
1031 // nucleon) is not taken into account.
1032 G4Nucleon* aNucleon = 0;
1033 nucleus->StartLoop();
1034 while ( ( aNucleon = nucleus->GetNextNucleon() ) ) { /* Loop checking, 07.08.2015, A.Ribon */
1035 nucleusMomentum += aNucleon->Get4Momentum();
1036 if ( aNucleon->AreYouHit() ) { // Involved nucleons
1037 // Consider in sumMasses the nominal, i.e. on-shell, masses of the nucleons
1038 // (not the current masses, which could be different because the nucleons are off-shell).
1039
1040 sumMasses += std::sqrt( sqr( aNucleon->GetDefinition()->GetPDGMass() )
1041 + aNucleon->Get4Momentum().perp2() );
1042 sumMasses += 20.0*MeV; // Separation energy for a nucleon
1043
1044 //residualExcitationEnergy += ExcitationEPerWoundedNucleon;
1045 residualExcitationEnergy += -ExcitationEPerWoundedNucleon*G4Log( G4UniformRand());
1046 residualMassNumber--;
1047 // The absolute value below is needed only in the case of anti-nucleus.
1048 residualCharge -= std::abs( G4int( aNucleon->GetDefinition()->GetPDGCharge() ) );
1049 } else { // Spectator nucleons
1050 residualMomentum += aNucleon->Get4Momentum();
1051 }
1052 }
1053 #ifdef debugPutOnMassShell
1054 G4cout << "ExcitationEnergyPerWoundedNucleon " << ExcitationEPerWoundedNucleon <<" MeV"<<G4endl
1055 << "\t Residual Charge, MassNumber " << residualCharge << " " << residualMassNumber
1056 << G4endl << "\t Initial Momentum " << nucleusMomentum/GeV<<" GeV"
1057 << G4endl << "\t Residual Momentum " << residualMomentum/GeV<<" GeV"<<G4endl;
1058 #endif
1059
1060 residualMomentum.setPz( 0.0 );
1061 residualMomentum.setE( 0.0 );
1062 if ( residualMassNumber == 0 ) {
1063 residualMass = 0.0;
1064 residualExcitationEnergy = 0.0;
1065 } else {
1067 GetIonMass( residualCharge, residualMassNumber );
1068 if ( residualMassNumber == 1 ) {
1069 residualExcitationEnergy = 0.0;
1070 }
1071 }
1072 sumMasses += std::sqrt( sqr( residualMass ) + residualMomentum.perp2() );
1073 return true;
1074}
1075
1076
1077//============================================================================
1078
1079G4bool G4QGSParticipants::
1080GenerateDeltaIsobar( const G4double sqrtS, // input parameter
1081 const G4int numberOfInvolvedNucleons, // input parameter
1082 G4Nucleon* involvedNucleons[], // input & output parameter
1083 G4double& sumMasses ) { // input & output parameter
1084
1085 // This method, which is called only by PutOnMassShell, check whether is possible to
1086 // re-interpret some of the involved nucleons as delta-isobars:
1087 // - either by replacing a proton (2212) with a Delta+ (2214),
1088 // - or by replacing a neutron (2112) with a Delta0 (2114).
1089 // The on-shell mass of these delta-isobars is ~1232 MeV, so ~292-294 MeV heavier than
1090 // the corresponding nucleon on-shell mass. However 400.0*MeV is considered to estimate
1091 // the max number of deltas compatible with the available energy.
1092 // The delta-isobars are considered with the same transverse momentum as their
1093 // corresponding nucleons.
1094 // This method assumes that all the parameters have been initialized by the caller;
1095 // the action of this method consists in modifying (eventually) involveNucleons and
1096 // sumMasses. The return value is "false" only in the case that the input parameters
1097 // have unphysical values.
1098
1099 if ( sqrtS < 0.0 || numberOfInvolvedNucleons <= 0 || sumMasses < 0.0 ) return false;
1100
1101 //const G4double ProbDeltaIsobar = 0.05; // Uzhi 6.07.2012
1102 //const G4double ProbDeltaIsobar = 0.25; // Uzhi 13.06.2013
1103 const G4double probDeltaIsobar = 0.10; // A.R. 07.08.2013
1104
1105 G4int maxNumberOfDeltas = G4int( (sqrtS - sumMasses)/(400.0*MeV) );
1106 G4int numberOfDeltas = 0;
1107
1108 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1109 //G4cout << "i maxNumberOfDeltas probDeltaIsobar " << i << " " << maxNumberOfDeltas
1110 // << " " << probDeltaIsobar << G4endl;
1111 if ( G4UniformRand() < probDeltaIsobar && numberOfDeltas < maxNumberOfDeltas ) {
1112 numberOfDeltas++;
1113 if ( ! involvedNucleons[i] ) continue;
1114 G4VSplitableHadron* splitableHadron = involvedNucleons[i]->GetSplitableHadron();
1115 G4double massNuc = std::sqrt( sqr( splitableHadron->GetDefinition()->GetPDGMass() )
1116 + splitableHadron->Get4Momentum().perp2() );
1117 //AR The absolute value below is needed in the case of an antinucleus.
1118 G4int pdgCode = std::abs( splitableHadron->GetDefinition()->GetPDGEncoding() );
1119 const G4ParticleDefinition* old_def = splitableHadron->GetDefinition();
1120 G4int newPdgCode = pdgCode/10; newPdgCode = newPdgCode*10 + 4; // Delta
1121 if ( splitableHadron->GetDefinition()->GetPDGEncoding() < 0 ) newPdgCode *= -1;
1122 const G4ParticleDefinition* ptr =
1124 splitableHadron->SetDefinition( ptr );
1125 G4double massDelta = std::sqrt( sqr( splitableHadron->GetDefinition()->GetPDGMass() )
1126 + splitableHadron->Get4Momentum().perp2() );
1127 //G4cout << i << " " << sqrtS/GeV << " " << sumMasses/GeV << " " << massDelta/GeV
1128 // << " " << massNuc << G4endl;
1129 if ( sqrtS < sumMasses + massDelta - massNuc ) { // Change cannot be accepted!
1130 splitableHadron->SetDefinition( old_def );
1131 break;
1132 } else { // Change is accepted
1133 sumMasses += ( massDelta - massNuc );
1134 }
1135 }
1136 }
1137 //G4cout << "maxNumberOfDeltas numberOfDeltas " << maxNumberOfDeltas << " "
1138 // << numberOfDeltas << G4endl;
1139 return true;
1140}
1141
1142
1143//============================================================================
1144
1145G4bool G4QGSParticipants::
1146SamplingNucleonKinematics( G4double averagePt2, // input parameter
1147 const G4double maxPt2, // input parameter
1148 G4double dCor, // input parameter
1149 G4V3DNucleus* nucleus, // input parameter
1150 const G4LorentzVector& pResidual, // input parameter
1151 const G4double residualMass, // input parameter
1152 const G4int residualMassNumber, // input parameter
1153 const G4int numberOfInvolvedNucleons, // input parameter
1154 G4Nucleon* involvedNucleons[], // input & output parameter
1155 G4double& mass2 ) { // output parameter
1156
1157 // This method, which is called only by PutOnMassShell, does the sampling of:
1158 // - either the target nucleons: this for any kind of hadronic interactions
1159 // (hadron-nucleus, nucleus-nucleus, antinucleus-nucleus);
1160 // - or the projectile nucleons or antinucleons: this only in the case of
1161 // nucleus-nucleus or antinucleus-nucleus interactions, respectively.
1162 // This method assumes that all the parameters have been initialized by the caller;
1163 // the action of this method consists in changing the properties of the nucleons
1164 // whose pointers are in the vector involvedNucleons, as well as changing the
1165 // variable mass2.
1166
1167 if ( ! nucleus ) return false;
1168
1169 if ( residualMassNumber == 0 && numberOfInvolvedNucleons == 1 ) {
1170 dCor = 0.0;
1171 averagePt2 = 0.0;
1172 }
1173
1174 G4bool success = true;
1175
1176 G4double SumMasses = residualMass;
1177 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1178 G4Nucleon* aNucleon = involvedNucleons[i];
1179 if ( ! aNucleon ) continue;
1180 SumMasses += aNucleon->GetSplitableHadron()->GetDefinition()->GetPDGMass();
1181 }
1182
1183 const G4int maxNumberOfLoops = 1000;
1184 G4int loopCounter = 0;
1185 do {
1186
1187 success = true;
1188 G4ThreeVector ptSum( 0.0, 0.0, 0.0 );
1189 G4double xSum = 0.0;
1190
1191 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1192 G4Nucleon* aNucleon = involvedNucleons[i];
1193 if ( ! aNucleon ) continue;
1194 G4ThreeVector tmpPt = GaussianPt( averagePt2, maxPt2 );
1195 ptSum += tmpPt;
1196 G4ThreeVector tmpX = GaussianPt( dCor*dCor, 1.0 );
1197 G4double x = tmpX.x() +
1198 aNucleon->GetSplitableHadron()->GetDefinition()->GetPDGMass()/SumMasses;
1199 if ( x < 0.0 || x > 1.0 ) {
1200 success = false;
1201 break;
1202 }
1203 xSum += x;
1204 //AR The energy is in the lab (instead of cms) frame but it will not be used.
1205 G4LorentzVector tmp( tmpPt.x(), tmpPt.y(), x, aNucleon->Get4Momentum().e() );
1206 aNucleon->SetMomentum( tmp );
1207 }
1208
1209 if ( xSum < 0.0 || xSum > 1.0 ) success = false;
1210
1211 if ( ! success ) continue;
1212
1213 G4double deltaPx = ( ptSum.x() - pResidual.x() ) / numberOfInvolvedNucleons;
1214 G4double deltaPy = ( ptSum.y() - pResidual.y() ) / numberOfInvolvedNucleons;
1215 G4double delta = 0.0;
1216 if ( residualMassNumber == 0 ) {
1217 delta = ( xSum - 1.0 ) / numberOfInvolvedNucleons;
1218 } else {
1219 delta = 0.0;
1220 }
1221
1222 xSum = 1.0;
1223 mass2 = 0.0;
1224 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1225 G4Nucleon* aNucleon = involvedNucleons[i];
1226 if ( ! aNucleon ) continue;
1227 G4double x = aNucleon->Get4Momentum().pz() - delta;
1228 xSum -= x;
1229 if ( residualMassNumber == 0 ) {
1230 if ( x <= 0.0 || x > 1.0 ) {
1231 success = false;
1232 break;
1233 }
1234 } else {
1235 if ( x <= 0.0 || x > 1.0 || xSum <= 0.0 || xSum > 1.0 ) {
1236 success = false;
1237 break;
1238 }
1239 }
1240 G4double px = aNucleon->Get4Momentum().px() - deltaPx;
1241 G4double py = aNucleon->Get4Momentum().py() - deltaPy;
1242 mass2 += ( sqr( aNucleon->GetSplitableHadron()->GetDefinition()->GetPDGMass() )
1243 + sqr( px ) + sqr( py ) ) / x;
1244 G4LorentzVector tmp( px, py, x, aNucleon->Get4Momentum().e() );
1245 aNucleon->SetMomentum( tmp );
1246 }
1247
1248 if ( success && residualMassNumber != 0 ) {
1249 mass2 += ( sqr( residualMass ) + pResidual.perp2() ) / xSum;
1250 }
1251
1252 #ifdef debugPutOnMassShell
1253 G4cout << "success " << success << G4endl << " Mt " << std::sqrt( mass2 )/GeV << G4endl;
1254 #endif
1255
1256 } while ( ( ! success ) &&
1257 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
1258 if ( loopCounter >= maxNumberOfLoops ) {
1259 success = false;
1260 }
1261
1262 return success;
1263}
1264
1265
1266//============================================================================
1267
1268G4bool G4QGSParticipants::
1269CheckKinematics( const G4double sValue, // input parameter
1270 const G4double sqrtS, // input parameter
1271 const G4double projectileMass2, // input parameter
1272 const G4double targetMass2, // input parameter
1273 const G4double nucleusY, // input parameter
1274 const G4bool isProjectileNucleus, // input parameter
1275 const G4int numberOfInvolvedNucleons, // input parameter
1276 G4Nucleon* involvedNucleons[], // input parameter
1277 G4double& targetWminus, // output parameter
1278 G4double& projectileWplus, // output parameter
1279 G4bool& success ) { // input & output parameter
1280
1281 // This method, which is called only by PutOnMassShell, checks whether the
1282 // kinematics is acceptable or not.
1283 // This method assumes that all the parameters have been initialized by the caller;
1284 // notice that the input boolean parameter isProjectileNucleus is meant to be true
1285 // only in the case of nucleus or antinucleus projectile.
1286 // The action of this method consists in computing targetWminus and projectileWplus
1287 // and setting the parameter success to false in the case that the kinematics should
1288 // be rejeted.
1289
1290 G4double decayMomentum2 = sqr( sValue ) + sqr( projectileMass2 ) + sqr( targetMass2 )
1291 - 2.0*sValue*projectileMass2 - 2.0*sValue*targetMass2
1292 - 2.0*projectileMass2*targetMass2;
1293 targetWminus = ( sValue - projectileMass2 + targetMass2 + std::sqrt( decayMomentum2 ) )
1294 / 2.0 / sqrtS;
1295 projectileWplus = sqrtS - targetMass2/targetWminus;
1296 G4double projectilePz = projectileWplus/2.0 - projectileMass2/2.0/projectileWplus;
1297 G4double projectileE = projectileWplus/2.0 + projectileMass2/2.0/projectileWplus;
1298 G4double projectileY(1.0e5);
1299 if (projectileE - projectilePz > 0.) {
1300 projectileY = 0.5 * G4Log( (projectileE + projectilePz)/
1301 (projectileE - projectilePz) );
1302 }
1303 G4double targetPz = -targetWminus/2.0 + targetMass2/2.0/targetWminus;
1304 G4double targetE = targetWminus/2.0 + targetMass2/2.0/targetWminus;
1305 G4double targetY = 0.5 * G4Log( (targetE + targetPz)/(targetE - targetPz) );
1306
1307 #ifdef debugPutOnMassShell
1308 G4cout << "decayMomentum2 " << decayMomentum2 << G4endl
1309 << "\t targetWminus projectileWplus " << targetWminus << " " << projectileWplus << G4endl
1310 << "\t projectileY targetY " << projectileY << " " << targetY << G4endl;
1311 #endif
1312
1313 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1314 G4Nucleon* aNucleon = involvedNucleons[i];
1315 if ( ! aNucleon ) continue;
1316 G4LorentzVector tmp = aNucleon->Get4Momentum();
1317 G4double mt2 = sqr( tmp.x() ) + sqr( tmp.y() ) +
1318 sqr( aNucleon->GetSplitableHadron()->GetDefinition()->GetPDGMass() );
1319 G4double x = tmp.z();
1320 G4double pz = -targetWminus*x/2.0 + mt2/(2.0*targetWminus*x);
1321 G4double e = targetWminus*x/2.0 + mt2/(2.0*targetWminus*x);
1322 if ( isProjectileNucleus ) {
1323 pz = projectileWplus*x/2.0 - mt2/(2.0*projectileWplus*x);
1324 e = projectileWplus*x/2.0 + mt2/(2.0*projectileWplus*x);
1325 }
1326 G4double nucleonY = 0.5 * G4Log( (e + pz)/(e - pz) );
1327
1328 #ifdef debugPutOnMassShell
1329 G4cout << "i nY pY nY-AY AY " << i << " " << nucleonY << " " << projectileY <<G4endl;
1330 #endif
1331
1332 if ( std::abs( nucleonY - nucleusY ) > 2 ||
1333 ( isProjectileNucleus && targetY > nucleonY ) ||
1334 ( ! isProjectileNucleus && projectileY < nucleonY ) ) {
1335 success = false;
1336 break;
1337 }
1338 }
1339 return true;
1340}
1341
1342
1343//============================================================================
1344
1345G4bool G4QGSParticipants::
1346FinalizeKinematics( const G4double w, // input parameter
1347 const G4bool isProjectileNucleus, // input parameter
1348 const G4LorentzRotation& boostFromCmsToLab, // input parameter
1349 const G4double residualMass, // input parameter
1350 const G4int residualMassNumber, // input parameter
1351 const G4int numberOfInvolvedNucleons, // input parameter
1352 G4Nucleon* involvedNucleons[], // input & output parameter
1353 G4LorentzVector& residual4Momentum ) { // output parameter
1354
1355 // This method, which is called only by PutOnMassShell, finalizes the kinematics:
1356 // this method is called when we are sure that the sampling of the kinematics is
1357 // acceptable.
1358 // This method assumes that all the parameters have been initialized by the caller;
1359 // notice that the input boolean parameter isProjectileNucleus is meant to be true
1360 // only in the case of nucleus or antinucleus projectile: this information is needed
1361 // because the sign of pz (in the center-of-mass frame) in this case is opposite
1362 // with respect to the case of a normal hadron projectile.
1363 // The action of this method consists in modifying the momenta of the nucleons
1364 // (in the lab frame) and computing the residual 4-momentum (in the center-of-mass
1365 // frame).
1366
1367 G4ThreeVector residual3Momentum( 0.0, 0.0, 1.0 );
1368
1369 for ( G4int i = 0; i < numberOfInvolvedNucleons; i++ ) {
1370 G4Nucleon* aNucleon = involvedNucleons[i];
1371 if ( ! aNucleon ) continue;
1372 G4LorentzVector tmp = aNucleon->Get4Momentum();
1373 residual3Momentum -= tmp.vect();
1374 G4double mt2 = sqr( tmp.x() ) + sqr( tmp.y() ) +
1375 sqr( aNucleon->GetSplitableHadron()->GetDefinition()->GetPDGMass() );
1376 G4double x = tmp.z();
1377 G4double pz = -w * x / 2.0 + mt2 / ( 2.0 * w * x );
1378 G4double e = w * x / 2.0 + mt2 / ( 2.0 * w * x );
1379 // Reverse the sign of pz in the case of nucleus or antinucleus projectile
1380 if ( isProjectileNucleus ) pz *= -1.0;
1381 tmp.setPz( pz );
1382 tmp.setE( e );
1383 tmp.transform( boostFromCmsToLab );
1384 aNucleon->SetMomentum( tmp );
1385 G4VSplitableHadron* splitableHadron = aNucleon->GetSplitableHadron();
1386 splitableHadron->Set4Momentum( tmp );
1387 #ifdef debugPutOnMassShell
1388 G4cout << "Target involved nucleon No, name, 4Mom "
1389 << i<<" "<<aNucleon->GetDefinition()->GetParticleName()<<" "<<tmp<< G4endl;
1390 #endif
1391 }
1392
1393 G4double residualMt2 = sqr( residualMass ) + sqr( residual3Momentum.x() )
1394 + sqr( residual3Momentum.y() );
1395
1396 #ifdef debugPutOnMassShell
1397 G4cout <<G4endl<< "w residual3Momentum.z() " << w << " " << residual3Momentum.z() << G4endl;
1398 #endif
1399
1400 G4double residualPz = 0.0;
1401 G4double residualE = 0.0;
1402 if ( residualMassNumber != 0 ) {
1403 residualPz = -w * residual3Momentum.z() / 2.0 +
1404 residualMt2 / ( 2.0 * w * residual3Momentum.z() );
1405 residualE = w * residual3Momentum.z() / 2.0 +
1406 residualMt2 / ( 2.0 * w * residual3Momentum.z() );
1407 // Reverse the sign of residualPz in the case of nucleus or antinucleus projectile
1408 if ( isProjectileNucleus ) residualPz *= -1.0;
1409 }
1410
1411 residual4Momentum.setPx( residual3Momentum.x() );
1412 residual4Momentum.setPy( residual3Momentum.y() );
1413 residual4Momentum.setPz( residualPz );
1414 residual4Momentum.setE( residualE );
1415
1416 return true;
1417}
1418
1419//======================================================
1421{
1422 #ifdef debugQGSParticipants
1423 G4cout<<G4endl<<"PerformDiffractiveCollisions()......"<<G4endl
1424 <<"theInteractions.size() "<<theInteractions.size()<<G4endl;
1425 #endif
1426
1427 unsigned int i;
1428 for (i = 0; i < theInteractions.size(); i++)
1429 {
1430 G4InteractionContent* anIniteraction = theInteractions[i];
1431 #ifdef debugQGSParticipants
1432 G4cout<<"Interaction # and its status "
1433 <<i<<" "<<theInteractions[i]->GetStatus()<<G4endl;
1434 #endif
1435
1436 G4int InterStatus = theInteractions[i]->GetStatus();
1437 if ( (InterStatus == PrD) || (InterStatus == TrD) || (InterStatus == DD))
1438 { // Selection of diffractive interactions
1439 #ifdef debugQGSParticipants
1440 G4cout<<"Simulation of diffractive interaction. "<<InterStatus <<" PrD/TrD/DD/ND/Qech - 0,1,2,3,4"<<G4endl;
1441 #endif
1442
1443 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1444
1445 #ifdef debugQGSParticipants
1446 G4cout<<"The proj. before inter "
1449 G4cout<<"The targ. before inter " <<aTarget->Get4Momentum()<<" "
1450 <<aTarget->Get4Momentum().mag()<<G4endl;
1451 #endif
1452
1453 if ( InterStatus == PrD )
1455
1456 if ( InterStatus == TrD )
1458
1459 if ( InterStatus == DD )
1461
1462 #ifdef debugQGSParticipants
1463 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1465 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1466 <<aTarget->Get4Momentum().mag()<<G4endl;
1467 #endif
1468 }
1469
1470 if ( InterStatus == Qexc )
1471 { // Quark exchange process
1472 #ifdef debugQGSParticipants
1473 G4cout<<"Simulation of interaction with quark exchange."<<G4endl;
1474 #endif
1475 G4VSplitableHadron* aTarget = anIniteraction->GetTarget();
1476
1477 #ifdef debugQGSParticipants
1478 G4cout<<"The proj. before inter " <<theProjectileSplitable->Get4Momentum()<<" "
1480 G4cout<<"The targ. before inter "<<aTarget->Get4Momentum()<<" "
1481 <<aTarget->Get4Momentum().mag()<<G4endl;
1482 #endif
1483
1485
1486 #ifdef debugQGSParticipants
1487 G4cout<<"The proj. after inter " <<theProjectileSplitable->Get4Momentum()<<" "
1489 G4cout<<"The targ. after inter " <<aTarget->Get4Momentum()<<" "
1490 <<aTarget->Get4Momentum().mag()<<G4endl;
1491 #endif
1492 }
1493 }
1494}
1495
1496//======================================================
1498{
1499 if ( ! theProjectileSplitable ) return false;
1500
1501 const G4double aHugeValue = 1.0e+10;
1502
1503 #ifdef debugQGSParticipants
1504 G4cout<<G4endl<<"DeterminePartonMomenta()......"<<G4endl;
1505 G4cout<<"theProjectile status (0 -nondiffr, #0 diffr./reggeon): "<<theProjectileSplitable->GetStatus()<<G4endl;
1506 #endif
1507
1508 if (theProjectileSplitable->GetStatus() != 0) {return false;} // There were only diffractive interactions.
1509
1510 G4LorentzVector Projectile4Momentum = theProjectileSplitable->Get4Momentum();
1511 G4LorentzVector Psum = Projectile4Momentum;
1512
1513 G4double VqM_pr(0.), VaqM_pr(0.), VqM_tr(350.), VqqM_tr(700);
1514 if (std::abs(theProjectile.GetDefinition()->GetBaryonNumber()) != 0) {VqM_pr=350*MeV; VaqM_pr=700*MeV;}
1515
1516 #ifdef debugQGSParticipants
1517 G4cout<<"Projectile 4 momentum "<<Psum<<G4endl
1518 <<"Target nucleon momenta at start"<<G4endl;
1519 #endif
1520
1521 std::vector<G4VSplitableHadron*>::iterator i;
1522 G4int NuclNo=0;
1523
1524 for (i = theTargets.begin(); i != theTargets.end(); i++ )
1525 {
1526 Psum += (*i)->Get4Momentum();
1527 #ifdef debugQGSParticipants
1528 G4cout<<"Nusleus nucleon # and its 4Mom. "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1529 #endif
1530 NuclNo++;
1531 }
1532
1533 G4LorentzRotation toCms( -1*Psum.boostVector() );
1534
1535 G4LorentzVector Ptmp = toCms*Projectile4Momentum;
1536
1537 toCms.rotateZ( -1*Ptmp.phi() );
1538 toCms.rotateY( -1*Ptmp.theta() );
1539 G4LorentzRotation toLab(toCms.inverse());
1540 Projectile4Momentum.transform( toCms );
1541 // Ptarget.transform( toCms );
1542
1543 #ifdef debugQGSParticipants
1544 G4cout<<G4endl<<"In CMS---------------"<<G4endl;
1545 G4cout<<"Projectile 4 Mom "<<Projectile4Momentum<<G4endl;
1546 #endif
1547
1548 NuclNo=0;
1549 G4LorentzVector Target4Momentum(0.,0.,0.,0.);
1550 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1551 {
1552 G4LorentzVector tmp= (*i)->Get4Momentum(); tmp.transform( toCms );
1553 (*i)->Set4Momentum( tmp );
1554 #ifdef debugQGSParticipants
1555 G4cout<<"Target nucleon # and 4Mom "<<" "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1556 #endif
1557 Target4Momentum += tmp;
1558 NuclNo++;
1559 }
1560
1561 G4double S = Psum.mag2();
1562 G4double SqrtS = std::sqrt(S);
1563
1564 #ifdef debugQGSParticipants
1565 G4cout<<"Sum of target nucleons 4 momentum "<<Target4Momentum<<G4endl<<G4endl;
1566 G4cout<<"Target nucleons mom: px, py, z_1, m_i"<<G4endl;
1567 #endif
1568
1569 //G4double PplusProjectile = Projectile4Momentum.plus();
1570 G4double PminusTarget = Target4Momentum.minus();
1571 NuclNo=0;
1572
1573 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1574 {
1575 G4LorentzVector tmp = (*i)->Get4Momentum(); // tmp.boost(bstToCM);
1576
1577 //AR-19Jan2017 : the following line is causing a strange crash when Geant4
1578 // is built in optimized mode.
1579 // To fix it, I get the mass square instead of directly the
1580 // mass from the Lorentz vector, and then I take care of the
1581 // square root. If the mass square is negative, a JustWarning
1582 // exception is thrown, and the mass is set to 0.
1583 //G4double Mass = tmp.mag();
1584 G4double Mass2 = tmp.mag2();
1585 G4double Mass = 0.0;
1586 if ( Mass2 < 0.0 ) {
1588 ed << "Projectile " << theProjectile.GetDefinition()->GetParticleName()
1589 << "  4-momentum " << Psum << G4endl;
1590 ed << "LorentzVector tmp " << tmp << " with Mass2 " << Mass2 << G4endl;
1591 G4Exception( "G4QGSParticipants::DeterminePartonMomenta(): 4-momentum with negative mass!",
1592 "HAD_QGSPARTICIPANTS_001", JustWarning, ed );
1593 } else {
1594 Mass = std::sqrt( Mass2 );
1595 }
1596
1597 tmp.setPz(tmp.minus()/PminusTarget); tmp.setE(Mass);
1598 (*i)->Set4Momentum(tmp);
1599 #ifdef debugQGSParticipants
1600 G4cout<<"Target nucleons # and mom: "<<NuclNo<<" "<<(*i)->Get4Momentum()<<G4endl;
1601 #endif
1602 NuclNo++;
1603 }
1604
1605 //+++++++++++++++++++++++++++++++++++++++++++
1606
1607 G4double SigPt = sigmaPt;
1608 G4Parton* aParton(0);
1609 G4ThreeVector aPtVector(0.,0.,0.);
1610 G4LorentzVector tmp(0.,0.,0.,0.);
1611
1612 G4double Mt(0.);
1613 G4double ProjSumMt(0.), ProjSumMt2perX(0.);
1614 G4double TargSumMt(0.), TargSumMt2perX(0.);
1615
1616
1617 G4double aBeta = beta; // Member of the class
1618
1619 const G4ParticleDefinition* theProjectileDefinition = theProjectileSplitable->GetDefinition();
1620 if (theProjectileDefinition == G4PionMinus::PionMinusDefinition()) aBeta = 1.;
1621 if (theProjectileDefinition == G4Gamma::GammaDefinition()) aBeta = 1.;
1622 if (theProjectileDefinition == G4PionPlus::PionPlusDefinition()) aBeta = 1.;
1623 if (theProjectileDefinition == G4PionZero::PionZeroDefinition()) aBeta = 1.;
1624 if (theProjectileDefinition == G4KaonPlus::KaonPlusDefinition()) aBeta = 0.;
1625 if (theProjectileDefinition == G4KaonMinus::KaonMinusDefinition()) aBeta = 0.;
1626
1627 G4double Xmin = 0.;
1628
1629 G4bool Success = true; G4int attempt = 0;
1630 const G4int maxNumberOfAttempts = 1000;
1631 do
1632 {
1633 attempt++; if( attempt == 100*(attempt/100) ) {SigPt/=2.;}
1634
1635 ProjSumMt=0.; ProjSumMt2perX=0.;
1636 TargSumMt=0.; TargSumMt2perX=0.;
1637
1638 Success = true;
1640 #ifdef debugQGSParticipants
1641 G4cout<<"attempt ------------------------ "<<attempt<<G4endl;
1642 G4cout<<"nSeaPair of proj "<<nSeaPair<<G4endl;
1643 #endif
1644
1645 G4double SumPx = 0.;
1646 G4double SumPy = 0.;
1647 G4double SumZ = 0.;
1648 G4int NumberOfUnsampledSeaQuarks = 2*nSeaPair;
1649
1650 G4double Qmass=0.;
1651 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1652 {
1653 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1654 #ifdef debugQGSParticipants
1655 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1656 #endif
1657 aPtVector = GaussianPt(SigPt, aHugeValue);
1658 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1659 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1660 Mt = std::sqrt(aPtVector.mag2()+sqr(Qmass));
1661 ProjSumMt += Mt;
1662
1663 // Sampling of Z fraction
1664 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1665 SumZ += tmp.z();
1666
1667 NumberOfUnsampledSeaQuarks--;
1668 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1669 tmp.setE(sqr(Mt));
1670 aParton->Set4Momentum(tmp);
1671
1672 aParton = theProjectileSplitable->GetNextAntiParton(); // for anti-quarks
1673 #ifdef debugQGSParticipants
1674 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1675 G4cout<<" "<<tmp<<" "<<SumZ<<G4endl;
1676 #endif
1677 aPtVector = GaussianPt(SigPt, aHugeValue);
1678 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1679 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1680 Mt = std::sqrt(aPtVector.mag2()+sqr(Qmass));
1681 ProjSumMt += Mt;
1682
1683 // Sampling of Z fraction
1684 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1685 SumZ += tmp.z();
1686
1687 NumberOfUnsampledSeaQuarks--;
1688 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1689 tmp.setE(sqr(Mt));
1690 aParton->Set4Momentum(tmp);
1691 #ifdef debugQGSParticipants
1692 G4cout<<" "<<tmp<<" "<<SumZ<<G4endl;
1693 #endif
1694 }
1695
1696 // For valence quark
1697 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1698 #ifdef debugQGSParticipants
1699 G4cout<<"Val quark of Pr"<<" "<<aParton->GetDefinition()->GetParticleName();
1700 #endif
1701 aPtVector = GaussianPt(SigPt, aHugeValue);
1702 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1703 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1704 Mt = std::sqrt(aPtVector.mag2()+sqr(VqM_pr));
1705 ProjSumMt += Mt;
1706
1707 // Sampling of Z fraction
1708 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1709 SumZ += tmp.z();
1710
1711 ProjSumMt2perX +=sqr(Mt)/tmp.pz();
1712 tmp.setE(sqr(Mt));
1713 aParton->Set4Momentum(tmp);
1714
1715 // For valence di-quark
1717 #ifdef debugQGSParticipants
1718 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1719 G4cout<<" "<<tmp<<" "<<SumZ<<" (z-fraction)"<<G4endl;
1720 #endif
1721 tmp.setPx(-SumPx); tmp.setPy(-SumPy);
1722 //Uzhi 2019 Mt = std::sqrt(aPtVector.mag2()+sqr(VaqM_pr));
1723 Mt = std::sqrt( sqr(SumPx) + sqr(SumPy) + sqr(VaqM_pr) ); //Uzhi 2019
1724 ProjSumMt += Mt;
1725 tmp.setPz(1.-SumZ);
1726
1727 ProjSumMt2perX +=sqr(Mt)/tmp.pz(); // QQmass=750 MeV
1728 tmp.setE(sqr(Mt));
1729 aParton->Set4Momentum(tmp);
1730 #ifdef debugQGSParticipants
1731 G4cout<<" "<<tmp<<" "<<SumZ+(1.-SumZ)<<" (z-fraction)"<<G4endl;
1732 #endif
1733
1734 // End of work with the projectile
1735
1736 // Work with target nucleons
1737
1738 NuclNo=0;
1739 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1740 {
1741 nSeaPair = (*i)->GetSoftCollisionCount()-1;
1742 #ifdef debugQGSParticipants
1743 G4cout<<"nSeaPair of target N "<<nSeaPair<<G4endl
1744 <<"Target nucleon 4Mom "<<(*i)->Get4Momentum()<<G4endl;
1745 #endif
1746
1747 SumPx = (*i)->Get4Momentum().px() * (-1.);
1748 SumPy = (*i)->Get4Momentum().py() * (-1.);
1749 SumZ = 0.;
1750
1751 G4double SumZw=0.;
1752 NumberOfUnsampledSeaQuarks = 2*nSeaPair;
1753
1754 Qmass=0;
1755 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1756 {
1757 aParton = (*i)->GetNextParton(); // for quarks
1758 #ifdef debugQGSParticipants
1759 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1760 #endif
1761 aPtVector = GaussianPt(SigPt, aHugeValue);
1762 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1763 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1764 Mt=std::sqrt(aPtVector.mag2()+sqr(Qmass));
1765 TargSumMt += Mt;
1766
1767 // Sampling of Z fraction
1768 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1769 SumZ += tmp.z();
1770 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1771 SumZw+=tmp.pz();
1772 NumberOfUnsampledSeaQuarks--;
1773 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1774 tmp.setE(sqr(Mt));
1775 aParton->Set4Momentum(tmp);
1776
1777 aParton = (*i)->GetNextAntiParton(); // for anti-quarks
1778 #ifdef debugQGSParticipants
1779 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1780 G4cout<<" "<<tmp<<" "<<SumZw<<" "<<SumZ<<G4endl;
1781 #endif
1782 aPtVector = GaussianPt(SigPt, aHugeValue);
1783 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1784 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1785 Mt=std::sqrt(aPtVector.mag2()+sqr(Qmass));
1786 TargSumMt += Mt;
1787
1788 // Sampling of Z fraction
1789 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1790 SumZ += tmp.z();
1791 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1792 SumZw+=tmp.pz();
1793 NumberOfUnsampledSeaQuarks--;
1794 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1795 tmp.setE(sqr(Mt));
1796 aParton->Set4Momentum(tmp);
1797 #ifdef debugQGSParticipants
1798 G4cout<<" "<<tmp<<" "<<SumZw<<" "<<SumZ<<G4endl;
1799 #endif
1800 }
1801
1802 // Valence quark
1803 aParton = (*i)->GetNextParton(); // for quarks
1804 #ifdef debugQGSParticipants
1805 G4cout<<"Val quark of Tr"<<" "<<aParton->GetDefinition()->GetParticleName();
1806 #endif
1807 aPtVector = GaussianPt(SigPt, aHugeValue);
1808 tmp.setPx(aPtVector.x()); tmp.setPy(aPtVector.y());
1809 SumPx += aPtVector.x(); SumPy += aPtVector.y();
1810 Mt=std::sqrt(aPtVector.mag2()+sqr(VqM_tr));
1811 TargSumMt += Mt;
1812
1813 // Sampling of Z fraction
1814 tmp.setPz(SampleX(Xmin, NumberOfUnsampledSeaQuarks, 2*nSeaPair, aBeta)*(1.0-SumZ));
1815 SumZ += tmp.z();
1816 tmp.setPz((*i)->Get4Momentum().pz()*tmp.pz());
1817 SumZw+=tmp.pz();
1818 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1819 tmp.setE(sqr(Mt));
1820 aParton->Set4Momentum(tmp);
1821
1822 // Valence di-quark
1823 aParton = (*i)->GetNextAntiParton(); // for quarks
1824 #ifdef debugQGSParticipants
1825 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1826 G4cout<<" "<<tmp<<" "<<SumZw<<" (sum z-fracs) "<<SumZ<<" (total z-sum) "<<G4endl;
1827 #endif
1828 tmp.setPx(-SumPx); tmp.setPy(-SumPy);
1829 //Uzhi 2019 Mt=std::sqrt(aPtVector.mag2()+sqr(VqqM_tr));
1830 Mt=std::sqrt( sqr(SumPx) + sqr(SumPy) + sqr(VqqM_tr) ); //Uzhi 2019
1831 TargSumMt += Mt;
1832
1833 tmp.setPz((*i)->Get4Momentum().pz()*(1.0 - SumZ));
1834 SumZw+=tmp.pz();
1835 TargSumMt2perX +=sqr(Mt)/tmp.pz();
1836 tmp.setE(sqr(Mt));
1837 aParton->Set4Momentum(tmp);
1838 #ifdef debugQGSParticipants
1839 G4cout<<" "<<tmp<<" "<<SumZw<<" "<<1.0<<" "<<(*i)->Get4Momentum().pz()<<G4endl;
1840 #endif
1841
1842 } // End of for(i = theTargets.begin(); i != theTargets.end(); i++ )
1843
1844 if( ProjSumMt + TargSumMt > SqrtS ) {
1845 Success = false; continue;}
1846 if( std::sqrt(ProjSumMt2perX) + std::sqrt(TargSumMt2perX) > SqrtS ) {
1847 Success = false; continue;}
1848
1849 } while( (!Success) &&
1850 attempt < maxNumberOfAttempts ); /* Loop checking, 07.08.2015, A.Ribon */
1851
1852 if ( attempt >= maxNumberOfAttempts ) {
1853 return false;
1854 }
1855
1856 //+++++++++++++++++++++++++++++++++++++++++++
1857
1858 G4double DecayMomentum2 = sqr(S) + sqr(ProjSumMt2perX) + sqr(TargSumMt2perX)
1859 - 2.0*S*ProjSumMt2perX - 2.0*S*TargSumMt2perX - 2.0*ProjSumMt2perX*TargSumMt2perX;
1860
1861 G4double targetWminus=( S - ProjSumMt2perX + TargSumMt2perX + std::sqrt( DecayMomentum2 ))/2.0/SqrtS;
1862 G4double projectileWplus = SqrtS - TargSumMt2perX/targetWminus;
1863
1864 G4LorentzVector Tmp(0.,0.,0.,0.);
1865 G4double z(0.);
1866
1868 #ifdef debugQGSParticipants
1869 G4cout<<"Backward transformation ===================="<<G4endl;
1870 G4cout<<"nSeaPair of proj "<<nSeaPair<<G4endl;
1871 #endif
1872
1873 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1874 {
1875 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1876 #ifdef debugQGSParticipants
1877 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1878 #endif
1879 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1880
1881 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1882 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1883 Tmp.transform( toLab );
1884
1885 aParton->Set4Momentum(Tmp);
1886
1887 aParton = theProjectileSplitable->GetNextAntiParton(); // for anti-quarks
1888 #ifdef debugQGSParticipants
1889 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1890 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1891 #endif
1892 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1893 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1894 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1895 Tmp.transform( toLab );
1896
1897 aParton->Set4Momentum(Tmp);
1898 #ifdef debugQGSParticipants
1899 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1900 #endif
1901 }
1902
1903 // For valence quark
1904 aParton = theProjectileSplitable->GetNextParton(); // for quarks
1905 #ifdef debugQGSParticipants
1906 G4cout<<"Val quark of Pr"<<" "<<aParton->GetDefinition()->GetParticleName();
1907 #endif
1908 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1909 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1910 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1911 Tmp.transform( toLab );
1912
1913 aParton->Set4Momentum(Tmp);
1914
1915 // For valence di-quark
1917 #ifdef debugQGSParticipants
1918 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1919 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1920 #endif
1921 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1922 Tmp.setPz(projectileWplus*z/2.0 - Tmp.e()/(2.0*z*projectileWplus));
1923 Tmp.setE( projectileWplus*z/2.0 + Tmp.e()/(2.0*z*projectileWplus));
1924 Tmp.transform( toLab );
1925
1926 aParton->Set4Momentum(Tmp);
1927
1928 #ifdef debugQGSParticipants
1929 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1930 #endif
1931
1932 // End of work with the projectile
1933
1934 // Work with target nucleons
1935 NuclNo=0;
1936 for(i = theTargets.begin(); i != theTargets.end(); i++ )
1937 {
1938 nSeaPair = (*i)->GetSoftCollisionCount()-1;
1939 #ifdef debugQGSParticipants
1940 G4cout<<"nSeaPair of target and N# "<<nSeaPair<<" "<<NuclNo<<G4endl;
1941 #endif
1942 NuclNo++;
1943 for (G4int aSeaPair = 0; aSeaPair < nSeaPair; aSeaPair++)
1944 {
1945 aParton = (*i)->GetNextParton(); // for quarks
1946 #ifdef debugQGSParticipants
1947 G4cout<<"Sea quarks: "<<aSeaPair<<" "<<aParton->GetDefinition()->GetParticleName();
1948 #endif
1949 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1950 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1951 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1952 Tmp.transform( toLab );
1953
1954 aParton->Set4Momentum(Tmp);
1955
1956 aParton = (*i)->GetNextAntiParton(); // for quarks
1957 #ifdef debugQGSParticipants
1958 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1959 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1960 #endif
1961 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1962 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1963 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1964 Tmp.transform( toLab );
1965
1966 aParton->Set4Momentum(Tmp);
1967 #ifdef debugQGSParticipants
1968 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<G4endl;
1969 #endif
1970 }
1971
1972 // Valence quark
1973
1974 aParton = (*i)->GetNextParton(); // for quarks
1975 #ifdef debugQGSParticipants
1976 G4cout<<"Val quark of Tr"<<" "<<aParton->GetDefinition()->GetParticleName();
1977 #endif
1978 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1979 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1980 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1981 Tmp.transform( toLab );
1982
1983 aParton->Set4Momentum(Tmp);
1984
1985 // Valence di-quark
1986 aParton = (*i)->GetNextAntiParton(); // for quarks
1987 #ifdef debugQGSParticipants
1988 G4cout<<" "<<aParton->GetDefinition()->GetParticleName()<<G4endl;
1989 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1990 #endif
1991 Tmp =aParton->Get4Momentum(); z=Tmp.z();
1992 Tmp.setPz(-targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1993 Tmp.setE( targetWminus*z/2.0 + Tmp.e()/(2.0*z*targetWminus));
1994 Tmp.transform( toLab );
1995
1996 aParton->Set4Momentum(Tmp);
1997 #ifdef debugQGSParticipants
1998 G4cout<<" "<<Tmp<<" "<<Tmp.mag()<<" (mass)"<<G4endl;
1999 #endif
2000 NuclNo++;
2001 } // End of for(i = theTargets.begin(); i != theTargets.end(); i++ )
2002
2003 return true;
2004}
2005
2006//======================================================
2008SampleX(G4double anXmin, G4int nSea, G4int totalSea, G4double aBeta)
2009{
2010 G4double Xmin=anXmin; G4int Nsea=totalSea; Xmin*=1.; Nsea++; // Must be erased
2011 G4double Oalfa = 1./(alpha + 1.);
2012 G4double Obeta = 1./(aBeta + (alpha + 1.)*nSea + 1.); // ?
2013
2014 G4double Ksi1, Ksi2, r1, r2, r12;
2015 const G4int maxNumberOfLoops = 1000;
2016 G4int loopCounter = 0;
2017 do
2018 {
2019 Ksi1 = G4UniformRand(); r1 = G4Pow::GetInstance()->powA(Ksi1,Oalfa);
2020 Ksi2 = G4UniformRand(); r2 = G4Pow::GetInstance()->powA(Ksi2,Obeta);
2021 r12=r1+r2;
2022 } while( ( r12 > 1.) &&
2023 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
2024 if ( loopCounter >= maxNumberOfLoops ) {
2025 return 0.5; // Just an acceptable value, without any physics consideration.
2026 }
2027
2028 G4double result = r1/r12;
2029 return result;
2030}
2031
2032//======================================================
2033void G4QGSParticipants::CreateStrings()
2034{
2035
2036 #ifdef debugQGSParticipants
2037 G4cout<<"CreateStrings() ..................."<<G4endl;
2038 #endif
2039
2040 if ( ! theProjectileSplitable ) {
2041 #ifdef debugQGSParticipants
2042 G4cout<<"BAD situation: theProjectileSplitable is NULL ! Returning immediately"<<G4endl;
2043 #endif
2044 return;
2045 }
2046
2047 #ifdef debugQGSParticipants
2048 G4cout<<"theProjectileSplitable->GetStatus() "<<theProjectileSplitable->GetStatus()<<G4endl;
2049 G4LorentzVector str4Mom;
2050 #endif
2051
2052 if( theProjectileSplitable->GetStatus() == 1 ) // The projectile has participated only in diffr. inter.
2053 {
2055
2059 #ifdef debugQGSParticipants
2060 G4cout << "Pr. Diffr. String: Qs 4mom X " <<G4endl;
2061 G4cout << " " << aPair->GetParton1()->GetPDGcode() << " "
2062 << aPair->GetParton1()->Get4Momentum() << " "
2063 << aPair->GetParton1()->GetX() << " " << G4endl;
2064 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2065 << aPair->GetParton2()->Get4Momentum() << " "
2066 << aPair->GetParton2()->GetX() << " " << G4endl;
2067 str4Mom += aPair->GetParton1()->Get4Momentum();
2068 str4Mom += aPair->GetParton2()->Get4Momentum();
2069 #endif
2070
2071 thePartonPairs.push_back(aPair);
2072 }
2073
2074 G4int N_EnvTarg = NumberOfInvolvedNucleonsOfTarget;
2075
2076 for ( G4int i = 0; i < N_EnvTarg; i++ ) {
2077 G4Nucleon* aTargetNucleon = TheInvolvedNucleonsOfTarget[ i ];
2078
2079 G4VSplitableHadron* HitTargetNucleon = aTargetNucleon->GetSplitableHadron();
2080
2081 #ifdef debugQGSParticipants
2082 G4cout<<"Involved Nucleon # and its status "<<i<<" "<<HitTargetNucleon->GetStatus()<<G4endl;
2083 #endif
2084 if( HitTargetNucleon->GetStatus() >= 1) // Create diffractive string
2085 {
2086 G4ThreeVector Position = HitTargetNucleon->GetPosition();
2087
2088 G4PartonPair * aPair = new G4PartonPair(HitTargetNucleon->GetNextParton(),
2089 HitTargetNucleon->GetNextAntiParton(),
2091 #ifdef debugQGSParticipants
2092 G4cout << "Tr. Diffr. String: Qs 4mom X " <<G4endl;
2093 G4cout << "Diffr. String " << aPair->GetParton1()->GetPDGcode() << " "
2094 << aPair->GetParton1()->Get4Momentum() << " "
2095 << aPair->GetParton1()->GetX() << " " << G4endl;
2096 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2097 << aPair->GetParton2()->Get4Momentum() << " "
2098 << aPair->GetParton2()->GetX() << " " << G4endl;
2099
2100 str4Mom += aPair->GetParton1()->Get4Momentum();
2101 str4Mom += aPair->GetParton2()->Get4Momentum();
2102 #endif
2103
2104 thePartonPairs.push_back(aPair);
2105 }
2106 }
2107
2108 //-----------------------------------------
2109 #ifdef debugQGSParticipants
2110 G4cout<<"Strings created in soft interactions"<<G4endl;
2111 #endif
2112 std::vector<G4InteractionContent*>::iterator i;
2113 G4int IntNo=0;
2114 i = theInteractions.begin();
2115 while ( i != theInteractions.end() ) /* Loop checking, 07.08.2015, A.Ribon */
2116 {
2117 G4InteractionContent* anIniteraction = *i;
2118 G4PartonPair * aPair = NULL;
2119
2120 #ifdef debugQGSParticipants
2121 G4cout<<"An interaction # and soft coll. # "<<IntNo<<" "
2122 <<anIniteraction->GetNumberOfSoftCollisions()<<G4endl;
2123 #endif
2124 IntNo++;
2125 if (anIniteraction->GetNumberOfSoftCollisions())
2126 {
2127 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
2128 G4VSplitableHadron* pTarget = anIniteraction->GetTarget();
2129
2130 for (G4int j = 0; j < anIniteraction->GetNumberOfSoftCollisions(); j++)
2131 {
2132 aPair = new G4PartonPair(pTarget->GetNextParton(), pProjectile->GetNextAntiParton(),
2134 #ifdef debugQGSParticipants
2135 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2136 << aPair->GetParton1()->Get4Momentum() << " "
2137 << aPair->GetParton1()->Get4Momentum().mag()<<G4endl;
2138 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2139 << aPair->GetParton2()->Get4Momentum() << " "
2140 <<aPair->GetParton2()->Get4Momentum().mag()<<G4endl;
2141 str4Mom += aPair->GetParton1()->Get4Momentum();
2142 str4Mom += aPair->GetParton2()->Get4Momentum();
2143 #endif
2144
2145 thePartonPairs.push_back(aPair);
2146
2147 aPair = new G4PartonPair(pProjectile->GetNextParton(), pTarget->GetNextAntiParton(),
2149 #ifdef debugQGSParticipants
2150 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2151 << aPair->GetParton1()->Get4Momentum() << " "
2152 << aPair->GetParton1()->Get4Momentum().mag()<<G4endl;
2153 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2154 << aPair->GetParton2()->Get4Momentum() << " "
2155 << aPair->GetParton2()->Get4Momentum().mag()<<G4endl;
2156 #endif
2157 #ifdef debugQGSParticipants
2158 str4Mom += aPair->GetParton1()->Get4Momentum();
2159 str4Mom += aPair->GetParton2()->Get4Momentum();
2160 #endif
2161
2162 thePartonPairs.push_back(aPair);
2163 }
2164
2165 delete *i;
2166 i=theInteractions.erase(i); // i now points to the next interaction
2167 }
2168 else
2169 {
2170 i++;
2171 }
2172 } // End of while ( i != theInteractions.end() )
2173 #ifdef debugQGSParticipants
2174 G4cout << "Sum of strings 4 momenta " << str4Mom << G4endl<<G4endl;
2175 #endif
2176}
2177
2178//============================================================================
2179
2180void G4QGSParticipants::GetResiduals() {
2181 // This method is needed for the correct application of G4PrecompoundModelInterface
2182
2183 #ifdef debugQGSParticipants
2184 G4cout << "GetResiduals(): GetProjectileNucleus()? " << GetProjectileNucleus() << G4endl;
2185 #endif
2186
2187 #ifdef debugQGSParticipants
2188 G4cout << "NumberOfInvolvedNucleonsOfTarget "<< NumberOfInvolvedNucleonsOfTarget << G4endl;
2189 #endif
2190
2191 G4double DeltaExcitationE = TargetResidualExcitationEnergy / G4double( NumberOfInvolvedNucleonsOfTarget );
2192 G4LorentzVector DeltaPResidualNucleus = TargetResidual4Momentum /
2193 G4double( NumberOfInvolvedNucleonsOfTarget );
2194
2195 for ( G4int i = 0; i < NumberOfInvolvedNucleonsOfTarget; i++ ) {
2196 G4Nucleon* aNucleon = TheInvolvedNucleonsOfTarget[i];
2197
2198 #ifdef debugQGSParticipants
2199 G4VSplitableHadron* targetSplitable = aNucleon->GetSplitableHadron();
2200 G4cout << i << " Hit? " << aNucleon->AreYouHit() << " " << targetSplitable << G4endl;
2201 if ( targetSplitable ) G4cout << i << "Status " << targetSplitable->GetStatus() << G4endl;
2202 #endif
2203
2204 G4LorentzVector tmp = -DeltaPResidualNucleus;
2205 aNucleon->SetMomentum( tmp );
2206 aNucleon->SetBindingEnergy( DeltaExcitationE );
2207 }
2208
2209 //-------------------------------------
2210 if( TargetResidualMassNumber != 0 )
2211 {
2212 G4ThreeVector bstToCM =TargetResidual4Momentum.findBoostToCM();
2213
2214 G4V3DNucleus* theTargetNucleus = GetTargetNucleus();
2215 G4LorentzVector residualMomentum(0.,0.,0.,0.);
2216 G4Nucleon* aNucleon = 0;
2217 theTargetNucleus->StartLoop();
2218 while ( ( aNucleon = theTargetNucleus->GetNextNucleon() ) ) { /* Loop checking, 07.08.2015, A.Ribon */
2219 if ( !aNucleon->AreYouHit() ) {
2220 G4LorentzVector tmp=aNucleon->Get4Momentum(); tmp.boost(bstToCM);
2221 aNucleon->SetMomentum(tmp);
2222 residualMomentum +=tmp;
2223 }
2224 }
2225
2226 residualMomentum/=TargetResidualMassNumber;
2227
2228 G4double Mass = TargetResidual4Momentum.mag();
2229 G4double SumMasses=0.;
2230
2231 aNucleon = 0;
2232 theTargetNucleus->StartLoop();
2233 while ( ( aNucleon = theTargetNucleus->GetNextNucleon() ) ) { /* Loop checking, 07.08.2015, A.Ribon */
2234 if ( !aNucleon->AreYouHit() ) {
2235 G4LorentzVector tmp=aNucleon->Get4Momentum() - residualMomentum;
2236 G4double E=std::sqrt(tmp.vect().mag2()+
2237 sqr(aNucleon->GetDefinition()->GetPDGMass()-aNucleon->GetBindingEnergy()));
2238 tmp.setE(E); aNucleon->SetMomentum(tmp);
2239 SumMasses+=E;
2240 }
2241 }
2242
2243 G4double Chigh=Mass/SumMasses; G4double Clow=0; G4double C;
2244 const G4int maxNumberOfLoops = 1000;
2245 G4int loopCounter = 0;
2246 do
2247 {
2248 C=(Chigh+Clow)/2.;
2249
2250 SumMasses=0.;
2251 aNucleon = 0;
2252 theTargetNucleus->StartLoop();
2253 while ( ( aNucleon = theTargetNucleus->GetNextNucleon() ) ) { /* Loop checking, 07.08.2015, A.Ribon */
2254 if ( !aNucleon->AreYouHit() ) {
2255 G4LorentzVector tmp=aNucleon->Get4Momentum();
2256 G4double E=std::sqrt(tmp.vect().mag2()*sqr(C)+
2257 sqr(aNucleon->GetDefinition()->GetPDGMass()-aNucleon->GetBindingEnergy()));
2258 SumMasses+=E;
2259 }
2260 }
2261
2262 if(SumMasses > Mass) {Chigh=C;}
2263 else {Clow =C;}
2264
2265 } while( (Chigh-Clow > 0.01) &&
2266 ++loopCounter < maxNumberOfLoops ); /* Loop checking, 07.08.2015, A.Ribon */
2267 if ( loopCounter >= maxNumberOfLoops ) {
2268 #ifdef debugQGSParticipants
2269 G4cout <<"BAD situation: forced loop exit!" << G4endl;
2270 #endif
2271 // Perhaps there is something to set here...
2272 } else {
2273 aNucleon = 0;
2274 theTargetNucleus->StartLoop();
2275 while ( ( aNucleon = theTargetNucleus->GetNextNucleon() ) ) { /* Loop checking, 07.08.2015, A.Ribon */
2276 if ( !aNucleon->AreYouHit() ) {
2277 G4LorentzVector tmp=aNucleon->Get4Momentum()*C;
2278 G4double E=std::sqrt(tmp.vect().mag2()+
2279 sqr(aNucleon->GetDefinition()->GetPDGMass()-aNucleon->GetBindingEnergy()));
2280 tmp.setE(E); tmp.boost(-bstToCM);
2281 aNucleon->SetMomentum(tmp);
2282 }
2283 }
2284 }
2285
2286 } // End of if( TargetResidualMassNumber != 0 )
2287 //-------------------------------------
2288
2289 #ifdef debugQGSParticipants
2290 G4cout << "End GetResiduals -----------------" << G4endl;
2291 #endif
2292
2293}
2294
2295
2296//======================================================
2298{
2299 std::vector<G4InteractionContent*>::iterator i;
2300 G4LorentzVector str4Mom;
2301 i = theInteractions.begin();
2302 while ( i != theInteractions.end() ) /* Loop checking, 07.08.2015, A.Ribon */
2303 {
2304 G4InteractionContent* anIniteraction = *i;
2305 G4PartonPair * aPair = NULL;
2306 if (anIniteraction->GetNumberOfSoftCollisions())
2307 {
2308 G4VSplitableHadron* pProjectile = anIniteraction->GetProjectile();
2309 G4VSplitableHadron* pTarget = anIniteraction->GetTarget();
2310 for (G4int j = 0; j < anIniteraction->GetNumberOfSoftCollisions(); j++)
2311 {
2312 aPair = new G4PartonPair(pTarget->GetNextParton(), pProjectile->GetNextAntiParton(),
2314 #ifdef debugQGSParticipants
2315 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2316 << aPair->GetParton1()->Get4Momentum() << " "
2317 << aPair->GetParton1()->GetX() << " " << G4endl;
2318 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2319 << aPair->GetParton2()->Get4Momentum() << " "
2320 << aPair->GetParton2()->GetX() << " " << G4endl;
2321 #endif
2322 #ifdef debugQGSParticipants
2323 str4Mom += aPair->GetParton1()->Get4Momentum();
2324 str4Mom += aPair->GetParton2()->Get4Momentum();
2325 #endif
2326 thePartonPairs.push_back(aPair);
2327 aPair = new G4PartonPair(pProjectile->GetNextParton(), pTarget->GetNextAntiParton(),
2329 #ifdef debugQGSParticipants
2330 G4cout << "SoftPair " << aPair->GetParton1()->GetPDGcode() << " "
2331 << aPair->GetParton1()->Get4Momentum() << " "
2332 << aPair->GetParton1()->GetX() << " " << G4endl;
2333 G4cout << " " << aPair->GetParton2()->GetPDGcode() << " "
2334 << aPair->GetParton2()->Get4Momentum() << " "
2335 << aPair->GetParton2()->GetX() << " " << G4endl;
2336 #endif
2337 #ifdef debugQGSParticipants
2338 str4Mom += aPair->GetParton1()->Get4Momentum();
2339 str4Mom += aPair->GetParton2()->Get4Momentum();
2340 #endif
2341 thePartonPairs.push_back(aPair);
2342 }
2343 delete *i;
2344 i=theInteractions.erase(i); // i now points to the next interaction
2345 } else {
2346 i++;
2347 }
2348 }
2349 #ifdef debugQGSParticipants
2350 G4cout << " string 4 mom " << str4Mom << G4endl;
2351 #endif
2352}
2353
2354
2355//************************************************
2357{
2358 // Check reaction threshold - goes to CheckThreshold
2359
2362
2363 G4LorentzVector aPrimaryMomentum(thePrimary.GetMomentum(), thePrimary.GetTotalEnergy());
2364 G4LorentzVector aTargetNMomentum(0.,0.,0.,938.);
2365
2366 if ((!(aPrimaryMomentum.e()>-1)) && (!(aPrimaryMomentum.e()<1)) )
2367 {
2368 throw G4HadronicException(__FILE__, __LINE__,
2369 "G4GammaParticipants::SelectInteractions: primary nan energy.");
2370 }
2371 G4double S = (aPrimaryMomentum + aTargetNMomentum).mag2();
2372 G4double ThresholdMass = thePrimary.GetMass() + 938.;
2373 ModelMode = SOFT;
2374
2375 #ifdef debugQGSParticipants
2376 G4cout << "Gamma projectile - SelectInteractions " << G4endl;
2377 G4cout<<"Energy and Nucleus "<<thePrimary.GetTotalEnergy()<<" "<<theNucleus->GetMassNumber()<<G4endl;
2378 G4cout << "SqrtS ThresholdMass ModelMode " <<std::sqrt(S)<<" "<<ThresholdMass<<" "<<ModelMode<< G4endl;
2379 #endif
2380
2381 if (sqr(ThresholdMass + ThresholdParameter) > S)
2382 {
2384 }
2385 if (sqr(ThresholdMass + QGSMThreshold) > S) // thus only diffractive in cascade!
2386 {
2388 }
2389
2390 // first find the collisions HPW
2391 std::for_each(theInteractions.begin(), theInteractions.end(), DeleteInteractionContent());
2392 theInteractions.clear();
2393 G4int totalCuts = 0;
2394
2396 G4int NucleonNo=0;
2397
2399 G4Nucleon * pNucleon = 0;
2400
2401 while( (pNucleon = theNucleus->GetNextNucleon()) ) /* Loop checking, 07.08.2015, A.Ribon */
2402 { if(NucleonNo == theCurrent) break; NucleonNo++; }
2403
2404 if ( pNucleon ) {
2405 G4QGSMSplitableHadron* aTarget = new G4QGSMSplitableHadron(*pNucleon);
2406
2407 pNucleon->Hit(aTarget);
2408
2409 if ( (0.06 > G4UniformRand() &&(ModelMode==SOFT)) || (ModelMode==DIFFRACTIVE ) )
2410 {
2413
2414 aInteraction->SetTarget(aTarget);
2415 aInteraction->SetTargetNucleon(pNucleon);
2416 aTarget->SetCollisionCount(0);
2417 aTarget->SetStatus(1);
2418
2419 aInteraction->SetNumberOfDiffractiveCollisions(1);
2420 aInteraction->SetNumberOfSoftCollisions(0);
2421 aInteraction->SetStatus(1);
2422
2423 theInteractions.push_back(aInteraction);
2424 totalCuts += 1;
2425 }
2426 else
2427 {
2428 // nondiffractive soft interaction occurs
2429 aTarget->IncrementCollisionCount(1);
2430 aTarget->SetStatus(0);
2431 theTargets.push_back(aTarget);
2432
2435
2437 aInteraction->SetTarget(aTarget);
2438 aInteraction->SetTargetNucleon(pNucleon);
2439 aInteraction->SetNumberOfSoftCollisions(1);
2440 aInteraction->SetStatus(3);
2441 theInteractions.push_back(aInteraction);
2442 totalCuts += 1;
2443 }
2444 }
2446}
double S(double temp)
double C(double temp)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ThreadLocal G4int G4QGSParticipants_NPart
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
#define G4UniformRand()
Definition: Randomize.hh:52
double x() const
double mag2() const
double y() const
double mag() const
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
double theta() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
double minus() const
double perp2() const
void setVect(const Hep3Vector &)
Hep3Vector findBoostToCM() const
HepLorentzVector & transform(const HepRotation &)
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:80
void SetNumberOfDiffractiveCollisions(int)
void SetTargetNucleon(G4Nucleon *aNucleon)
G4Nucleon * GetTargetNucleon() const
G4VSplitableHadron * GetProjectile() const
void SetTarget(G4VSplitableHadron *aTarget)
void SetStatus(G4int aValue)
G4VSplitableHadron * GetTarget() const
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:107
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:107
const G4ThreeVector & GetPosition() const
Definition: G4Nucleon.hh:138
G4bool AreYouHit() const
Definition: G4Nucleon.hh:96
void SetPosition(const G4ThreeVector aPosition)
Definition: G4Nucleon.hh:133
void SetMomentum(G4LorentzVector &aMomentum)
Definition: G4Nucleon.hh:70
G4VSplitableHadron * GetSplitableHadron() const
Definition: G4Nucleon.hh:95
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
void Hit(G4VSplitableHadron *aHit)
Definition: G4Nucleon.hh:89
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
void SetBindingEnergy(G4double anEnergy)
Definition: G4Nucleon.hh:74
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:84
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4Parton * GetParton2(void)
Definition: G4PartonPair.hh:76
G4Parton * GetParton1(void)
Definition: G4PartonPair.hh:71
G4double GetX()
Definition: G4Parton.hh:87
G4ParticleDefinition * GetDefinition()
Definition: G4Parton.hh:161
G4int GetPDGcode() const
Definition: G4Parton.hh:127
void Set4Momentum(const G4LorentzVector &aMomentum)
Definition: G4Parton.hh:148
const G4LorentzVector & Get4Momentum() const
Definition: G4Parton.hh:143
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:92
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:92
static G4PionZero * PionZeroDefinition()
Definition: G4PionZero.cc:102
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4bool ProjectileDiffraction=TRUE) const
virtual G4Parton * GetNextParton()
virtual G4Parton * GetNextAntiParton()
const G4double ThresholdParameter
G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta)
G4ThreeVector theCurrentVelocity
G4QuarkExchange theQuarkExchange
G4SingleDiffractiveExcitation theSingleDiffExcitation
virtual void DoLorentzBoost(G4ThreeVector aBoost)
std::vector< G4InteractionContent * > theInteractions
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
G4QGSDiffractiveExcitation theDiffExcitaton
G4QGSMSplitableHadron * theProjectileSplitable
virtual ~G4QGSParticipants()
const G4double QGSMThreshold
const G4double theNucleonRadius
void BuildInteractions(const G4ReactionProduct &thePrimary)
std::vector< G4PartonPair * > thePartonPairs
std::vector< G4VSplitableHadron * > theTargets
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
G4double GetMass() const
G4int ncPomerons()
Definition: G4Reggeons.cc:525
void SetS(G4double S)
Definition: G4Reggeons.cc:341
void GetProbabilities(G4double B, G4int Mode, G4double &Pint, G4double &Pprd, G4double &Ptrd, G4double &Pdd, G4double &Pnd, G4double &Pnvr)
Definition: G4Reggeons.cc:448
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4bool ProjectileDiffraction) const
virtual void SortNucleonsIncZ()=0
virtual G4double GetOuterRadius()=0
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4bool StartLoop()=0
virtual void Init(G4int theA, G4int theZ)=0
std::pair< G4double, G4double > ChooseImpactXandY(G4double maxImpact)
Definition: G4V3DNucleus.hh:86
virtual G4int GetMassNumber()=0
virtual void SetProjectileNucleus(G4V3DNucleus *aNucleus)
G4V3DNucleus * theNucleus
void SetTimeOfCreation(G4double aTime)
void SetStatus(const G4int aStatus)
void SetCollisionCount(G4int aCount)
const G4ParticleDefinition * GetDefinition() const
void Set4Momentum(const G4LorentzVector &a4Momentum)
virtual G4Parton * GetNextParton()=0
virtual G4Parton * GetNextAntiParton()=0
void SetDefinition(const G4ParticleDefinition *aDefinition)
const G4LorentzVector & Get4Momentum() const
const G4ThreeVector & GetPosition() const
void IncrementCollisionCount(G4int aCount)
T sqr(const T &x)
Definition: templates.hh:128
#define G4ThreadLocal
Definition: tls.hh:77