Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/parton_string/diffraction/src/G4FTFParticipants.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/hadronic/models/parton_string/diffraction/src/G4FTFParticipants.cc (Version 11.3.0) and /processes/hadronic/models/parton_string/diffraction/src/G4FTFParticipants.cc (Version 10.3.p2)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
 26 //                                                 26 //
                                                   >>  27 // $Id: G4FTFParticipants.cc 100828 2016-11-02 15:25:59Z gcosmo $
                                                   >>  28 // GEANT4 tag $Name:  $
 27 //                                                 29 //
 28                                                    30 
 29 // -------------------------------------------     31 // ------------------------------------------------------------
 30 //      GEANT 4 class implementation file          32 //      GEANT 4 class implementation file
 31 //                                                 33 //
 32 //      ---------------- G4FTFParticipants----     34 //      ---------------- G4FTFParticipants----------------
 33 //             by Gunter Folger, June 1998.        35 //             by Gunter Folger, June 1998.
 34 //       class finding colliding particles in      36 //       class finding colliding particles in FTFPartonStringModel
 35 //  Changed in a part by V. Uzhinsky in oder t     37 //  Changed in a part by V. Uzhinsky in oder to put in correcpondence
 36 //        with original FRITIOF mode. November     38 //        with original FRITIOF mode. November - December 2006.
 37 //  Ajusted for (anti) nucleus - nucleus inter     39 //  Ajusted for (anti) nucleus - nucleus interactions by V. Uzhinsky.
 38 //                    (February 2011)              40 //                    (February 2011)
 39 // -------------------------------------------     41 // ------------------------------------------------------------
 40                                                    42 
 41 #include <utility>                                 43 #include <utility>
 42 #include <vector>                                  44 #include <vector>
 43 #include <algorithm>                               45 #include <algorithm>
 44                                                    46 
 45 #include "G4FTFParticipants.hh"                    47 #include "G4FTFParticipants.hh"
 46 #include "G4ios.hh"                                48 #include "G4ios.hh"
 47 #include "Randomize.hh"                            49 #include "Randomize.hh"
 48 #include "G4SystemOfUnits.hh"                      50 #include "G4SystemOfUnits.hh"
 49 #include "G4FTFParameters.hh"                      51 #include "G4FTFParameters.hh"                            
 50 #include "G4DiffractiveSplitableHadron.hh"         52 #include "G4DiffractiveSplitableHadron.hh"
 51 #include "G4VSplitableHadron.hh"                   53 #include "G4VSplitableHadron.hh"
 52 #include "G4PhysicalConstants.hh"              << 
 53                                                    54 
 54                                                    55 
 55 //============================================     56 //============================================================================
 56                                                    57 
 57 //#define debugFTFparticipant                      58 //#define debugFTFparticipant
 58                                                    59 
 59                                                    60 
 60 //============================================     61 //============================================================================
 61                                                    62 
 62 G4FTFParticipants::G4FTFParticipants() : Bimpa <<  63 G4FTFParticipants::G4FTFParticipants() : currentInteraction( -1 ) {}
 63                                          Bmin2 <<  64 
 64                                          curre <<  65 
 65 {}                                             <<  66 //============================================================================
                                                   >>  67 
                                                   >>  68 G4FTFParticipants::G4FTFParticipants( const G4FTFParticipants& ) : 
                                                   >>  69   G4VParticipants(), currentInteraction( -1 ) 
                                                   >>  70 {
                                                   >>  71   G4Exception( "G4FTFParticipants::G4FTFParticipants()", "HAD_FTF_001",
                                                   >>  72                FatalException, " Must not use copy ctor()" );
                                                   >>  73 }
 66                                                    74 
 67                                                    75 
 68 //============================================     76 //============================================================================
 69                                                    77 
 70 G4FTFParticipants::~G4FTFParticipants() {}         78 G4FTFParticipants::~G4FTFParticipants() {}
 71                                                    79 
                                                   >>  80 
 72 //============================================     81 //============================================================================
 73                                                    82 
 74 void G4FTFParticipants::GetList( const G4React     83 void G4FTFParticipants::GetList( const G4ReactionProduct& thePrimary, 
 75                                  G4FTFParamete     84                                  G4FTFParameters* theParameters ) { 
 76                                                    85 
 77   #ifdef debugFTFparticipant                       86   #ifdef debugFTFparticipant
 78   G4cout << "Participants::GetList" << G4endl      87   G4cout << "Participants::GetList" << G4endl 
 79          << "thePrimary " << thePrimary.GetMom     88          << "thePrimary " << thePrimary.GetMomentum() << G4endl << G4endl;
 80   #endif                                           89   #endif
 81                                                    90 
 82   G4double betta_z = thePrimary.GetMomentum().     91   G4double betta_z = thePrimary.GetMomentum().z() / thePrimary.GetTotalEnergy();
 83   if ( betta_z < 1.0e-10 ) betta_z = 1.0e-10;      92   if ( betta_z < 1.0e-10 ) betta_z = 1.0e-10;
 84                                                    93 
 85   StartLoop(); // reset Loop over Interactions     94   StartLoop(); // reset Loop over Interactions
 86                                                    95 
 87   for ( unsigned int i = 0; i < theInteraction     96   for ( unsigned int i = 0; i < theInteractions.size(); i++ ) delete theInteractions[i];
 88   theInteractions.clear();                         97   theInteractions.clear();
 89                                                    98 
 90   G4double deltaxy = 2.0 * fermi;  // Extra nu     99   G4double deltaxy = 2.0 * fermi;  // Extra nuclear radius
 91                                                   100 
 92   if ( theProjectileNucleus == nullptr ) {  // << 101   if ( theProjectileNucleus == 0 ) {  // Hadron-nucleus or anti-baryon-nucleus interactions
 93                                                   102 
 94     G4double impactX( 0.0 ), impactY( 0.0 );      103     G4double impactX( 0.0 ), impactY( 0.0 );
 95     G4double B( 0.0 ), B2( 0.0 );              << 
 96                                                   104 
 97     G4VSplitableHadron* primarySplitable = new    105     G4VSplitableHadron* primarySplitable = new G4DiffractiveSplitableHadron( thePrimary );
 98                                                   106 
 99     #ifdef debugFTFparticipant                    107     #ifdef debugFTFparticipant
100     G4cout << "Hadron-nucleus or anti-baryon-n    108     G4cout << "Hadron-nucleus or anti-baryon-nucleus interactions" << G4endl;
101     #endif                                        109     #endif
102                                                   110 
103     G4double xyradius = theNucleus->GetOuterRa << 111     G4double xyradius;                          
                                                   >> 112     xyradius = theNucleus->GetOuterRadius() + deltaxy; // Range of impact parameter sampling
104                                                   113 
105     const G4int maxNumberOfLoops = 1000;          114     const G4int maxNumberOfLoops = 1000;
106     G4int loopCounter = 0;                        115     G4int loopCounter = 0;                                          
107     do {                                          116     do {
108                                                   117 
109       std::pair< G4double, G4double > theImpac    118       std::pair< G4double, G4double > theImpactParameter;
110       if ( SampleBinInterval() ) {             << 119       theImpactParameter = theNucleus->ChooseImpactXandY( xyradius );
111         B2 = GetBmin2() + G4UniformRand() * (  << 120       impactX = theImpactParameter.first; 
112         B  = B2 > 0.0 ? std::sqrt( B2 ) : 0.0; << 121       impactY = theImpactParameter.second;
113         G4double Phi = twopi * G4UniformRand() << 
114         impactX = B * std::cos( Phi );         << 
115         impactY = B * std::sin( Phi );         << 
116         SetImpactParameter( B );               << 
117       } else {                                 << 
118         theImpactParameter = theNucleus->Choos << 
119         impactX = theImpactParameter.first;    << 
120         impactY = theImpactParameter.second;   << 
121         SetImpactParameter( std::sqrt( sqr(imp << 
122       }                                        << 
123                                                   122 
124       #ifdef debugFTFparticipant                  123       #ifdef debugFTFparticipant
125       G4cout << "New interaction list," << " b << 124       G4cout << "New interaction list," << " b= " 
126              << std::sqrt( sqr(impactX ) + sqr    125              << std::sqrt( sqr(impactX ) + sqr( impactY ) )/fermi << G4endl;
127       #endif                                      126       #endif
128                                                   127 
129       G4ThreeVector thePosition( impactX, impa    128       G4ThreeVector thePosition( impactX, impactY, 0.0 );
130       primarySplitable->SetPosition( thePositi    129       primarySplitable->SetPosition( thePosition );
131                                                   130 
132       theNucleus->StartLoop();                    131       theNucleus->StartLoop();
133       G4Nucleon* nucleon;                         132       G4Nucleon* nucleon;
134                                                   133 
135       #ifdef debugFTFparticipant                  134       #ifdef debugFTFparticipant
136       G4int TrN( 0 );                             135       G4int TrN( 0 );
137       #endif                                      136       #endif
138                                                   137 
139       while ( ( nucleon = theNucleus->GetNextN    138       while ( ( nucleon = theNucleus->GetNextNucleon() ) ) {  /* Loop checking, 10.08.2015, A.Ribon */
140                                                   139 
141         G4double impact2 = sqr( impactX - nucl    140         G4double impact2 = sqr( impactX - nucleon->GetPosition().x() ) +
142                            sqr( impactY - nucl    141                            sqr( impactY - nucleon->GetPosition().y() );
143                                                   142 
144         if ( theParameters->GetProbabilityOfIn    143         if ( theParameters->GetProbabilityOfInteraction( impact2/fermi/fermi ) > 
145              G4UniformRand() ) {                  144              G4UniformRand() ) {
146           primarySplitable->SetStatus( 1 );  /    145           primarySplitable->SetStatus( 1 );  // It takes part in the interaction
147           G4VSplitableHadron* targetSplitable     146           G4VSplitableHadron* targetSplitable = 0;
148           if ( ! nucleon->AreYouHit() ) {         147           if ( ! nucleon->AreYouHit() ) {
149             targetSplitable = new G4Diffractiv    148             targetSplitable = new G4DiffractiveSplitableHadron( *nucleon );
150             nucleon->Hit( targetSplitable );      149             nucleon->Hit( targetSplitable );
151             targetSplitable->SetStatus( 1 ); /    150             targetSplitable->SetStatus( 1 ); // It takes part in the interaction
152                                                   151 
153             #ifdef debugFTFparticipant            152             #ifdef debugFTFparticipant
154             G4cout << "Participated nucleons #    153             G4cout << "Participated nucleons #, " << TrN << " " << "Splitable Pr* Tr* "
155                    << primarySplitable << " "     154                    << primarySplitable << " " << targetSplitable << G4endl;
156             #endif                                155             #endif
157                                                   156 
158           }                                       157           }
159           G4InteractionContent* aInteraction =    158           G4InteractionContent* aInteraction = new G4InteractionContent( primarySplitable );
160           G4Nucleon* PrNucleon = 0;               159           G4Nucleon* PrNucleon = 0;
161           aInteraction->SetProjectileNucleon(     160           aInteraction->SetProjectileNucleon( PrNucleon );
162           aInteraction->SetTarget( targetSplit    161           aInteraction->SetTarget( targetSplitable );
163           aInteraction->SetTargetNucleon( nucl    162           aInteraction->SetTargetNucleon( nucleon );     
164           aInteraction->SetStatus( 1 );           163           aInteraction->SetStatus( 1 );                  
165           aInteraction->SetInteractionTime( (     164           aInteraction->SetInteractionTime( ( primarySplitable->GetPosition().z() + 
166                                                   165                                               nucleon->GetPosition().z() ) / betta_z );
167           theInteractions.push_back( aInteract    166           theInteractions.push_back( aInteraction );
168         }                                         167         }
169                                                   168 
170         #ifdef debugFTFparticipant                169         #ifdef debugFTFparticipant
171         TrN++;                                    170         TrN++;
172         #endif                                    171         #endif
173                                                   172 
174       }                                           173       } 
175                                                   174 
176     } while ( ( theInteractions.size() == 0 )     175     } while ( ( theInteractions.size() == 0 ) && 
177               ++loopCounter < maxNumberOfLoops    176               ++loopCounter < maxNumberOfLoops );  /* Loop checking, 10.08.2015, A.Ribon */
178     if ( loopCounter >= maxNumberOfLoops ) {      177     if ( loopCounter >= maxNumberOfLoops ) {
179       #ifdef debugFTFparticipant                  178       #ifdef debugFTFparticipant
180       G4cout << "BAD situation: forced exit fr    179       G4cout << "BAD situation: forced exit from the while loop!" << G4endl;
181       #endif                                      180       #endif
182       return;                                     181       return;
183     }                                             182     }
184                                                   183 
185     #ifdef debugFTFparticipant                    184     #ifdef debugFTFparticipant
186     G4cout << "Number of Hit nucleons " << the << 185     G4cout << "Number of Hit nucleons " << theInteractions.size() << "\t Bx " << impactX/fermi
187            << "\t By[fm] " << impactY/fermi << << 186            << "\t By " << impactY/fermi << "\t B " 
188            << std::sqrt( sqr( impactX ) + sqr(    187            << std::sqrt( sqr( impactX ) + sqr( impactY ) )/fermi << G4endl << G4endl;
189     #endif                                        188     #endif
190                                                   189 
191     //SortInteractionsIncT(); // Not need beca    190     //SortInteractionsIncT(); // Not need because nucleons are sorted in increasing z-coordinates.
192     ShiftInteractionTime();  // To put correct    191     ShiftInteractionTime();  // To put correct times and z-coordinates
193     return;                                       192     return;
194                                                   193 
195   }  // end of if ( theProjectileNucleus == 0     194   }  // end of if ( theProjectileNucleus == 0 )
196                                                   195 
197   // Projectile and target are nuclei             196   // Projectile and target are nuclei
198                                                   197 
199   #ifdef debugFTFparticipant                      198   #ifdef debugFTFparticipant
200   G4cout << "Projectile and target are nuclei"    199   G4cout << "Projectile and target are nuclei" << G4endl;
201   #endif                                          200   #endif
202                                                   201 
203   //G4cout<<theProjectileNucleus->GetOuterRadi << 202 //G4cout<<theProjectileNucleus->GetOuterRadius()/fermi<<" "<<theNucleus->GetOuterRadius()/fermi<<" "<<deltaxy/fermi<<G4endl;
204                                                   203 
205   // Range of impact parameter sampling        << 204   G4double xyradius;                          
206   G4double xyradius = theProjectileNucleus->Ge << 205   xyradius = theProjectileNucleus->GetOuterRadius() +  // Range of impact parameter sampling
                                                   >> 206              theNucleus->GetOuterRadius() + deltaxy;
207                                                   207 
208   G4double impactX( 0.0 ), impactY( 0.0 );        208   G4double impactX( 0.0 ), impactY( 0.0 );
209   G4double B( 0.0 ), B2( 0.0 );                << 
210                                                   209 
211   const G4int maxNumberOfLoops = 1000;            210   const G4int maxNumberOfLoops = 1000;
212   G4int loopCounter = 0;                          211   G4int loopCounter = 0;
213   do {                                            212   do {
214                                                   213 
215     std::pair< G4double, G4double > theImpactP    214     std::pair< G4double, G4double > theImpactParameter;
216     if ( SampleBinInterval() ) {               << 215     theImpactParameter = theNucleus->ChooseImpactXandY( xyradius );
217       B2 = GetBmin2() + G4UniformRand() * ( Ge << 216     impactX = theImpactParameter.first; 
218       B  = B2 > 0.0 ? std::sqrt( B2 ) : 0.0;   << 217     impactY = theImpactParameter.second;
219       G4double Phi = twopi * G4UniformRand();  << 
220       impactX = B * std::cos( Phi );           << 
221       impactY = B * std::sin( Phi );           << 
222       SetImpactParameter( B );                 << 
223     } else {                                   << 
224       theImpactParameter = theNucleus->ChooseI << 
225       impactX = theImpactParameter.first;      << 
226       impactY = theImpactParameter.second;     << 
227       SetImpactParameter( std::sqrt( sqr(impac << 
228     }                                          << 
229                                                   218 
230     #ifdef debugFTFparticipant                    219     #ifdef debugFTFparticipant
231     G4cout << "New interaction list, " << "b[f << 220     G4cout << "New interaction list, " << "b "
232            << std::sqrt( sqr( impactX ) + sqr(    221            << std::sqrt( sqr( impactX ) + sqr( impactY ) )/fermi << G4endl;
233     #endif                                        222     #endif
234                                                   223 
235     G4ThreeVector theBeamPosition( impactX, im    224     G4ThreeVector theBeamPosition( impactX, impactY, 0.0 );                    
236                                                   225 
237     theProjectileNucleus->StartLoop();            226     theProjectileNucleus->StartLoop();
238     G4Nucleon* ProjectileNucleon;                 227     G4Nucleon* ProjectileNucleon;
239                                                   228 
240     #ifdef debugFTFparticipant                    229     #ifdef debugFTFparticipant
241     G4int PrNuclN( 0 );                           230     G4int PrNuclN( 0 );
242     #endif                                        231     #endif
243                                                   232 
244     while ( ( ProjectileNucleon = theProjectil    233     while ( ( ProjectileNucleon = theProjectileNucleus->GetNextNucleon() ) ) {  /* Loop checking, 10.08.2015, A.Ribon */
245                                                   234  
246       G4VSplitableHadron* ProjectileSplitable     235       G4VSplitableHadron* ProjectileSplitable = 0;
247       theNucleus->StartLoop();                    236       theNucleus->StartLoop();
248       G4Nucleon* TargetNucleon = 0;               237       G4Nucleon* TargetNucleon = 0;
249                                                   238 
250       #ifdef debugFTFparticipant                  239       #ifdef debugFTFparticipant
251       G4int TrNuclN( 0 );                         240       G4int TrNuclN( 0 );
252       #endif                                      241       #endif
253                                                   242 
254       while ( ( TargetNucleon = theNucleus->Ge    243       while ( ( TargetNucleon = theNucleus->GetNextNucleon() ) ) {  /* Loop checking, 10.08.2015, A.Ribon */
255                                                   244 
256         G4double impact2 = sqr( impactX + Proj    245         G4double impact2 = sqr( impactX + ProjectileNucleon->GetPosition().x() - 
257                                 TargetNucleon-    246                                 TargetNucleon->GetPosition().x() ) +
258                            sqr( impactY + Proj    247                            sqr( impactY + ProjectileNucleon->GetPosition().y() -
259                                 TargetNucleon-    248                                 TargetNucleon->GetPosition().y() );
260         G4VSplitableHadron* TargetSplitable =     249         G4VSplitableHadron* TargetSplitable = 0;
261         if ( theParameters->GetProbabilityOfIn    250         if ( theParameters->GetProbabilityOfInteraction( impact2/fermi/fermi ) >
262              G4UniformRand() ) {  // An intera    251              G4UniformRand() ) {  // An interaction has happend!
263                                                   252 
264           #ifdef debugFTFparticipant              253           #ifdef debugFTFparticipant
265           G4cout << G4endl << "An Interaction     254           G4cout << G4endl << "An Interaction has happend" << G4endl << "Proj N mom " << PrNuclN
266                  << " " << ProjectileNucleon->    255                  << " " << ProjectileNucleon->Get4Momentum() << "-------------" << G4endl
267                  << "Targ N mom " << TrNuclN <    256                  << "Targ N mom " << TrNuclN << " " << TargetNucleon->Get4Momentum() << G4endl
268                  << "PrN TrN Z coords [fm]" << << 257                  << "PrN TrN Z coords " << ProjectileNucleon->GetPosition().z()/fermi 
269                  << " " << TargetNucleon->GetP    258                  << " " << TargetNucleon->GetPosition().z()/fermi 
270                  << " " << ProjectileNucleon->    259                  << " " << ProjectileNucleon->GetPosition().z()/fermi + 
271                            TargetNucleon->GetP    260                            TargetNucleon->GetPosition().z()/fermi << G4endl;
272           #endif                                  261           #endif
273                                                   262 
274           if ( ! ProjectileNucleon->AreYouHit(    263           if ( ! ProjectileNucleon->AreYouHit() ) { 
275             // Projectile nucleon was not invo    264             // Projectile nucleon was not involved until now.
276             ProjectileSplitable = new G4Diffra    265             ProjectileSplitable = new G4DiffractiveSplitableHadron( *ProjectileNucleon );
277             ProjectileNucleon->Hit( Projectile    266             ProjectileNucleon->Hit( ProjectileSplitable );
278             ProjectileSplitable->SetStatus( 1     267             ProjectileSplitable->SetStatus( 1 );  // It takes part in the interaction
279           } else {  // Projectile nucleon was     268           } else {  // Projectile nucleon was involved before.
280             ProjectileSplitable = ProjectileNu    269             ProjectileSplitable = ProjectileNucleon->GetSplitableHadron();
281           }                                       270           }
282                                                   271 
283           if ( ! TargetNucleon->AreYouHit() )     272           if ( ! TargetNucleon->AreYouHit() ) {  // Target nucleon was not involved until now
284             TargetSplitable = new G4Diffractiv    273             TargetSplitable = new G4DiffractiveSplitableHadron( *TargetNucleon );
285             TargetNucleon->Hit( TargetSplitabl    274             TargetNucleon->Hit( TargetSplitable );
286             TargetSplitable->SetStatus( 1 );      275             TargetSplitable->SetStatus( 1 );   // It takes part in the interaction
287           } else {  // Target nucleon was invo    276           } else {  // Target nucleon was involved before.
288             TargetSplitable = TargetNucleon->G    277             TargetSplitable = TargetNucleon->GetSplitableHadron();
289           }                                       278           }
290                                                   279 
291           G4InteractionContent* anInteraction     280           G4InteractionContent* anInteraction = new G4InteractionContent( ProjectileSplitable );
292           anInteraction->SetTarget( TargetSpli    281           anInteraction->SetTarget( TargetSplitable );
293           anInteraction->SetProjectileNucleon(    282           anInteraction->SetProjectileNucleon( ProjectileNucleon );
294           anInteraction->SetTargetNucleon( Tar    283           anInteraction->SetTargetNucleon( TargetNucleon );
295           anInteraction->SetInteractionTime( (    284           anInteraction->SetInteractionTime( ( ProjectileNucleon->GetPosition().z() + 
296                                                   285                                                TargetNucleon->GetPosition().z() ) / betta_z );
297           anInteraction->SetStatus( 1 );          286           anInteraction->SetStatus( 1 );                     
298                                                   287 
299           #ifdef debugFTFparticipant              288           #ifdef debugFTFparticipant
300           G4cout << "Part anInteraction->GetIn << 289           G4cout << "Part anInteraction->GetInteractionTime() " 
301                  << anInteraction->GetInteract    290                  << anInteraction->GetInteractionTime()/fermi << G4endl
302                  << "Splitable Pr* Tr* " << Pr    291                  << "Splitable Pr* Tr* " << ProjectileSplitable << " " 
303                  << TargetSplitable << G4endl;    292                  << TargetSplitable << G4endl;
304           #endif                                  293           #endif
305                                                   294 
306           theInteractions.push_back( anInterac    295           theInteractions.push_back( anInteraction );
307                                                   296 
308         } // End of an Interaction has happend    297         } // End of an Interaction has happend!
309                                                   298 
310         #ifdef debugFTFparticipant                299         #ifdef debugFTFparticipant
311         TrNuclN++;                                300         TrNuclN++;
312         #endif                                    301         #endif
313                                                   302 
314       }  // End of while ( ( TargetNucleon = t    303       }  // End of while ( ( TargetNucleon = theNucleus->GetNextNucleon() ) )
315                                                   304 
316       #ifdef debugFTFparticipant                  305       #ifdef debugFTFparticipant
317       PrNuclN++;                                  306       PrNuclN++;
318       #endif                                      307       #endif
319                                                   308 
320     }  // End of while ( ( ProjectileNucleon =    309     }  // End of while ( ( ProjectileNucleon = theProjectileNucleus->GetNextNucleon() ) )
321                                                   310 
322     if ( theInteractions.size() != 0 ) theProj    311     if ( theInteractions.size() != 0 ) theProjectileNucleus->DoTranslation( theBeamPosition );
323                                                   312 
324   } while ( ( theInteractions.size() == 0 ) &&    313   } while ( ( theInteractions.size() == 0 ) &&
325             ++loopCounter < maxNumberOfLoops )    314             ++loopCounter < maxNumberOfLoops );  /* Loop checking, 10.08.2015, A.Ribon */
326   if ( loopCounter >= maxNumberOfLoops ) {        315   if ( loopCounter >= maxNumberOfLoops ) {
327     #ifdef debugFTFparticipant                    316     #ifdef debugFTFparticipant
328     G4cout << "BAD situation: forced exit from    317     G4cout << "BAD situation: forced exit from the while loop!" << G4endl;
329     #endif                                        318     #endif
330     return;                                       319     return;
331   }                                               320   }
332                                                   321 
333   SortInteractionsIncT();                         322   SortInteractionsIncT();
334   ShiftInteractionTime();                         323   ShiftInteractionTime();
335                                                   324 
336   #ifdef debugFTFparticipant                      325   #ifdef debugFTFparticipant
337   G4cout << G4endl << "Number of primary colli    326   G4cout << G4endl << "Number of primary collisions " << theInteractions.size() 
338          << "\t Bx[fm] " << impactX/fermi << " << 327          << "\t Bx " << impactX/fermi << "\t By " << impactY/fermi
339          << "\t B[fm] " << std::sqrt( sqr( imp << 328          << "\t B " << std::sqrt( sqr( impactX ) + sqr( impactY ) )/fermi << G4endl
340          << "FTF participant End. ############    329          << "FTF participant End. #######################" << G4endl << G4endl;
341   #endif                                          330   #endif
342   return;                                         331   return;
343 }                                                 332 }
344                                                   333 
345                                                   334 
346 //============================================    335 //============================================================================
347                                                   336 
348 bool G4FTFPartHelperForSortInT( const G4Intera    337 bool G4FTFPartHelperForSortInT( const G4InteractionContent* Int1, 
349                                 const G4Intera    338                                 const G4InteractionContent* Int2 ) {
350   return Int1->GetInteractionTime() < Int2->Ge    339   return Int1->GetInteractionTime() < Int2->GetInteractionTime();
351 }                                                 340 }
352                                                   341 
353                                                   342 
354 //============================================    343 //============================================================================
355                                                   344 
356 void G4FTFParticipants::SortInteractionsIncT()    345 void G4FTFParticipants::SortInteractionsIncT() {  // on increased T 
357   if ( theInteractions.size() < 2 ) return;  /    346   if ( theInteractions.size() < 2 ) return;  // Avoid unnecesary work
358   std::sort( theInteractions.begin(), theInter    347   std::sort( theInteractions.begin(), theInteractions.end(), G4FTFPartHelperForSortInT ); 
359 }                                                 348 }
360                                                   349 
361                                                   350 
362 //============================================    351 //============================================================================
363                                                   352 
364 void G4FTFParticipants::ShiftInteractionTime()    353 void G4FTFParticipants::ShiftInteractionTime() {
365   G4double InitialTime = theInteractions[0]->G    354   G4double InitialTime = theInteractions[0]->GetInteractionTime();
366   for ( unsigned int i = 1; i < theInteraction    355   for ( unsigned int i = 1; i < theInteractions.size(); i++ ) {
367     G4double InterTime = theInteractions[i]->G    356     G4double InterTime = theInteractions[i]->GetInteractionTime() - InitialTime;
368     theInteractions[i]->SetInteractionTime( In    357     theInteractions[i]->SetInteractionTime( InterTime );
369     G4InteractionContent* aCollision = theInte    358     G4InteractionContent* aCollision = theInteractions[i];
370     G4VSplitableHadron* projectile = aCollisio    359     G4VSplitableHadron* projectile = aCollision->GetProjectile();
371     G4VSplitableHadron* target = aCollision->G    360     G4VSplitableHadron* target = aCollision->GetTarget();
372     G4ThreeVector prPosition = projectile->Get    361     G4ThreeVector prPosition = projectile->GetPosition();
373     prPosition.setZ( target->GetPosition().z()    362     prPosition.setZ( target->GetPosition().z() );
374     projectile->SetPosition( prPosition );        363     projectile->SetPosition( prPosition );
375     projectile->SetTimeOfCreation( InterTime )    364     projectile->SetTimeOfCreation( InterTime );
376     target->SetTimeOfCreation( InterTime );       365     target->SetTimeOfCreation( InterTime );
377   }                                               366   }
378   return;                                         367   return;
379 }                                                 368 }
380                                                   369 
381                                                   370 
382 //============================================    371 //============================================================================
383                                                   372 
384 void G4FTFParticipants::Clean() {                 373 void G4FTFParticipants::Clean() {
385   for ( size_t i = 0; i < theInteractions.size    374   for ( size_t i = 0; i < theInteractions.size(); i++ ) {
386     if ( theInteractions[ i ] ) {                 375     if ( theInteractions[ i ] ) {
387       delete theInteractions[ i ];                376       delete theInteractions[ i ];
388       theInteractions[ i ] = 0;                   377       theInteractions[ i ] = 0;
389     }                                             378     }
390   }                                               379   }
391   theInteractions.clear();                        380   theInteractions.clear();
392   currentInteraction = -1;                        381   currentInteraction = -1;
393 }                                                 382 }
394                                                   383 
395                                                   384