Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GeneratorPrecompoundInterface.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//
27// -----------------------------------------------------------------------------
28// GEANT 4 class file
29//
30// History: first implementation
31// HPW, 10DEC 98, the decay part originally written by Gunter Folger
32// in his FTF-test-program.
33//
34// M.Kelsey, 28 Jul 2011 -- Replace loop to decay input secondaries
35// with new utility class, simplify cleanup loops
36//
37// A.Ribon, 27 Oct 2021 -- Extended the method PropagateNuclNucl
38// to deal with projectile hypernuclei and anti-hypernuclei
39//
40// -----------------------------------------------------------------------------
41
42#include <algorithm>
43#include <vector>
44
47#include "G4SystemOfUnits.hh"
50#include "G4Proton.hh"
51#include "G4Neutron.hh"
52
53#include "G4Deuteron.hh"
54#include "G4Triton.hh"
55#include "G4He3.hh"
56#include "G4Alpha.hh"
57
58#include "G4V3DNucleus.hh"
59#include "G4Nucleon.hh"
60
61#include "G4AntiProton.hh"
62#include "G4AntiNeutron.hh"
63#include "G4AntiDeuteron.hh"
64#include "G4AntiTriton.hh"
65#include "G4AntiHe3.hh"
66#include "G4AntiAlpha.hh"
67
68#include "G4HyperTriton.hh"
69#include "G4HyperH4.hh"
70#include "G4HyperAlpha.hh"
71#include "G4HyperHe5.hh"
72#include "G4DoubleHyperH4.hh"
74
75#include "G4AntiHyperTriton.hh"
76#include "G4AntiHyperH4.hh"
77#include "G4AntiHyperAlpha.hh"
78#include "G4AntiHyperHe5.hh"
81
82#include "G4FragmentVector.hh"
83#include "G4ReactionProduct.hh"
85#include "G4PreCompoundModel.hh"
89
92//---------------------------------------------------------------------
93#include "Randomize.hh"
94#include "G4Log.hh"
95
96//#define debugPrecoInt
97
99 : CaptureThreshold(70*MeV), DeltaM(5.0*MeV), DeltaR(0.0), secID(-1)
100{
103
106 He3 =G4He3::He3();
108
111
116
117 if(preModel) { SetDeExcitation(preModel); }
118 else {
121 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(hadi);
122 if(!pre) { pre = new G4PreCompoundModel(); }
123 SetDeExcitation(pre);
124 }
125
127}
128
130{
131}
132
134Propagate(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus)
135{
136 #ifdef debugPrecoInt
137 G4cout<<G4endl<<"G4GeneratorPrecompoundInterface::Propagate"<<G4endl;
138 G4cout<<"Target A and Z "<<theNucleus->GetMassNumber()<<" "<<theNucleus->GetCharge()<<G4endl;
139 G4cout<<"Directly produced particles number "<<theSecondaries->size()<<G4endl;
140 #endif
141
143
144 // decay the strong resonances
145 G4DecayKineticTracks decay(theSecondaries);
146 #ifdef debugPrecoInt
147 G4cout<<"Final stable particles number "<<theSecondaries->size()<<G4endl;
148 #endif
149
150 // prepare the fragment (it is assumed that target nuclei are never hypernuclei)
151 G4int anA=theNucleus->GetMassNumber();
152 G4int aZ=theNucleus->GetCharge();
153// G4double TargetNucleusMass = G4NucleiProperties::GetNuclearMass(anA, aZ);
154
155 G4int numberOfEx = 0;
156 G4int numberOfCh = 0;
157 G4int numberOfHoles = 0;
158
159 G4double R = theNucleus->GetNuclearRadius();
160
161 G4LorentzVector captured4Momentum(0.,0.,0.,0.);
162 G4LorentzVector Residual4Momentum(0.,0.,0.,0.); // TargetNucleusMass is not need at the moment
163 G4LorentzVector Secondary4Momentum(0.,0.,0.,0.);
164
165 // loop over secondaries
166 G4KineticTrackVector::iterator iter;
167 for(iter=theSecondaries->begin(); iter !=theSecondaries->end(); ++iter)
168 {
169 const G4ParticleDefinition* part = (*iter)->GetDefinition();
170 G4double e = (*iter)->Get4Momentum().e();
171 G4double mass = (*iter)->Get4Momentum().mag();
172 G4ThreeVector mom = (*iter)->Get4Momentum().vect();
173 if((part != proton && part != neutron) ||
174 ((*iter)->GetPosition().mag() > R)) {
175 G4ReactionProduct * theNew = new G4ReactionProduct(part);
176 theNew->SetMomentum(mom);
177 theNew->SetTotalEnergy(e);
178 theNew->SetCreatorModelID((*iter)->GetCreatorModelID());
179 theTotalResult->push_back(theNew);
180 Secondary4Momentum += (*iter)->Get4Momentum();
181 #ifdef debugPrecoInt
182 G4cout<<"Secondary 4Mom "<<part->GetParticleName()<<" "<<(*iter)->Get4Momentum()<<" "
183 <<(*iter)->Get4Momentum().mag()<<G4endl;
184 #endif
185 } else {
186 if( e-mass > -CaptureThreshold*G4Log( G4UniformRand()) ) {
187 G4ReactionProduct * theNew = new G4ReactionProduct(part);
188 theNew->SetMomentum(mom);
189 theNew->SetTotalEnergy(e);
190 theNew->SetCreatorModelID((*iter)->GetCreatorModelID());
191 theTotalResult->push_back(theNew);
192 Secondary4Momentum += (*iter)->Get4Momentum();
193 #ifdef debugPrecoInt
194 G4cout<<"Secondary 4Mom "<<part->GetParticleName()<<" "<<(*iter)->Get4Momentum()<<" "
195 <<(*iter)->Get4Momentum().mag()<<G4endl;
196 #endif
197 } else {
198 // within the nucleus, neutron or proton
199 // now calculate A, Z of the fragment, momentum, number of exciton states
200 ++anA;
201 ++numberOfEx;
202 G4int Z = G4int(part->GetPDGCharge()/eplus + 0.1);
203 aZ += Z;
204 numberOfCh += Z;
205 captured4Momentum += (*iter)->Get4Momentum();
206 #ifdef debugPrecoInt
207 G4cout<<"Captured 4Mom "<<part->GetParticleName()<<(*iter)->Get4Momentum()<<G4endl;
208 #endif
209 }
210 }
211 delete (*iter);
212 }
213 delete theSecondaries;
214
215 // loop over wounded nucleus
216 G4Nucleon * theCurrentNucleon =
217 theNucleus->StartLoop() ? theNucleus->GetNextNucleon() : 0;
218 while(theCurrentNucleon) /* Loop checking, 31.08.2015, G.Folger */
219 {
220 if(theCurrentNucleon->AreYouHit()) {
221 ++numberOfHoles;
222 ++numberOfEx;
223 --anA;
224 aZ -= G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/eplus + 0.1);
225
226 Residual4Momentum -= theCurrentNucleon->Get4Momentum();
227 }
228 theCurrentNucleon = theNucleus->GetNextNucleon();
229 }
230
231 #ifdef debugPrecoInt
232 G4cout<<G4endl;
233 G4cout<<"Secondary 4Mom "<<Secondary4Momentum<<G4endl;
234 G4cout<<"Captured 4Mom "<<captured4Momentum<<G4endl;
235 G4cout<<"Sec + Captured "<<Secondary4Momentum+captured4Momentum<<G4endl;
236 G4cout<<"Residual4Mom "<<Residual4Momentum<<G4endl;
237 G4cout<<"Sum 4 momenta "
238 <<Secondary4Momentum + captured4Momentum + Residual4Momentum <<G4endl;
239 #endif
240
241 // Check that we use QGS model; loop over wounded nucleus
242 G4bool QGSM(false);
243 theCurrentNucleon = theNucleus->StartLoop() ? theNucleus->GetNextNucleon() : 0;
244 while(theCurrentNucleon) /* Loop checking, 31.08.2015, G.Folger */
245 {
246 if(theCurrentNucleon->AreYouHit())
247 {
248 if(theCurrentNucleon->Get4Momentum().mag() <
249 theCurrentNucleon->GetDefinition()->GetPDGMass()) QGSM=true;
250 }
251 theCurrentNucleon = theNucleus->GetNextNucleon();
252 }
253
254 #ifdef debugPrecoInt
255 if(!QGSM){
256 G4cout<<G4endl;
257 G4cout<<"Residual A and Z "<<anA<<" "<<aZ<<G4endl;
258 G4cout<<"Residual 4Mom "<<Residual4Momentum<<G4endl;
259 if(numberOfEx == 0)
260 {G4cout<<"Residual 4Mom = 0 means that there were not wounded and captured nucleons"<<G4endl;}
261 }
262 #endif
263
264 if(anA == 0) return theTotalResult;
265
266 G4LorentzVector exciton4Momentum(0.,0.,0.,0.);
267 if(anA >= aZ)
268 {
269 if(!QGSM)
270 { // FTF model was used
272
273// G4LorentzVector exciton4Momentum = Residual4Momentum + captured4Momentum;
274 exciton4Momentum = Residual4Momentum + captured4Momentum;
275//exciton4Momentum.setE(std::sqrt(exciton4Momentum.vect().mag2()+sqr(fMass)));
276 G4double ActualMass = exciton4Momentum.mag();
277 if(ActualMass <= fMass ) {
278 exciton4Momentum.setE(std::sqrt(exciton4Momentum.vect().mag2()+sqr(fMass)));
279 }
280
281 #ifdef debugPrecoInt
282 G4double exEnergy = 0.0;
283 if(ActualMass <= fMass ) {exEnergy = 0.;}
284 else {exEnergy = ActualMass - fMass;}
285 G4cout<<"Ground state residual Mass "<<fMass<<" E* "<<exEnergy<<G4endl;
286 #endif
287 }
288 else
289 { // QGS model was used
290 G4double InitialTargetMass =
291 G4NucleiProperties::GetNuclearMass(theNucleus->GetMassNumber(), theNucleus->GetCharge());
292
293 exciton4Momentum =
294 GetPrimaryProjectile()->Get4Momentum() + G4LorentzVector(0.,0.,0.,InitialTargetMass)
295 -Secondary4Momentum;
296
298 G4double ActualMass = exciton4Momentum.mag();
299
300 #ifdef debugPrecoInt
301 G4cout<<G4endl;
302 G4cout<<"Residual A and Z "<<anA<<" "<<aZ<<G4endl;
303 G4cout<<"Residual4Momentum "<<exciton4Momentum<<G4endl;
304 G4cout<<"ResidualMass, GroundStateMass and E* "<<ActualMass<<" "<<fMass<<" "
305 <<ActualMass - fMass<<G4endl;
306 #endif
307
308 if(ActualMass - fMass < 0.)
309 {
310 G4double ResE = std::sqrt(exciton4Momentum.vect().mag2() + sqr(fMass+10*MeV));
311 exciton4Momentum.setE(ResE);
312 #ifdef debugPrecoInt
313 G4cout<<"ActualMass - fMass < 0. "<<ActualMass<<" "<<fMass<<" "<<ActualMass - fMass<<G4endl;
314 G4int Uzhi; G4cin>>Uzhi;
315 #endif
316 }
317 }
318
319 // Need to de-excite the remnant nucleus only if excitation energy > 0.
320 G4Fragment anInitialState(anA, aZ, exciton4Momentum);
321 anInitialState.SetNumberOfParticles(numberOfEx-numberOfHoles);
322 anInitialState.SetNumberOfCharged(numberOfCh);
323 anInitialState.SetNumberOfHoles(numberOfHoles);
324 anInitialState.SetCreatorModelID(secID);
325
326 G4ReactionProductVector * aPrecoResult =
327 theDeExcitation->DeExcite(anInitialState);
328 // fill pre-compound part into the result, and return
329 #ifdef debugPrecoInt
330 G4cout<<"Target fragment number "<<aPrecoResult->size()<<G4endl;
331 #endif
332 for(unsigned int ll=0; ll<aPrecoResult->size(); ++ll)
333 {
334 theTotalResult->push_back(aPrecoResult->operator[](ll));
335 #ifdef debugPrecoInt
336 G4cout<<"Fragment "<<ll<<" "
337 <<aPrecoResult->operator[](ll)->GetDefinition()->GetParticleName()<<" "
338 <<aPrecoResult->operator[](ll)->GetMomentum()<<" "
339 <<aPrecoResult->operator[](ll)->GetTotalEnergy()<<" "
340 <<aPrecoResult->operator[](ll)->GetDefinition()->GetPDGMass()<<G4endl;
341 #endif
342 }
343 delete aPrecoResult;
344 }
345
346 return theTotalResult;
347}
348
351{
352 G4cout << "G4GeneratorPrecompoundInterface: ApplyYourself interface called stand-allone."
353 << G4endl;
354 G4cout << "This class is only a mediator between generator and precompound"<<G4endl;
355 G4cout << "Please remove from your physics list."<<G4endl;
356 throw G4HadronicException(__FILE__, __LINE__, "SEVERE: G4GeneratorPrecompoundInterface model interface called stand-allone.");
357 return new G4HadFinalState;
358}
360{
361 outFile << "G4GeneratorPrecompoundInterface interfaces a high\n"
362 << "energy model through the wounded nucleus to precompound de-excitation.\n"
363 << "Low energy protons and neutron present among secondaries produced by \n"
364 << "the high energy generator and within the nucleus are captured. The wounded\n"
365 << "nucleus and the captured particles form an excited nuclear fragment. This\n"
366 << "fragment is passed to the Geant4 pre-compound model for de-excitation.\n"
367 << "Nuclear de-excitation:\n";
368 // preco
369
370}
371
372
374PropagateNuclNucl(G4KineticTrackVector* theSecondaries, G4V3DNucleus* theNucleus,
375 G4V3DNucleus* theProjectileNucleus)
376{
377#ifdef debugPrecoInt
378 G4cout<<G4endl<<"G4GeneratorPrecompoundInterface::PropagateNuclNucl "<<G4endl;
379 G4cout<<"Projectile A and Z (and numberOfLambdas) "<<theProjectileNucleus->GetMassNumber()<<" "
380 <<theProjectileNucleus->GetCharge()<<" ("
381 <<theProjectileNucleus->GetNumberOfLambdas()<<")"<<G4endl;
382 G4cout<<"Target A and Z "<<theNucleus->GetMassNumber()<<" "
383 <<theNucleus->GetCharge()<<" ("
384 <<theNucleus->GetNumberOfLambdas()<<")"<<G4endl;
385 G4cout<<"Directly produced particles number "<<theSecondaries->size()<<G4endl;
386 G4cout<<"Projectile 4Mom and mass "<<GetPrimaryProjectile()->Get4Momentum()<<" "
388#endif
389
390 // prepare the target residual (assumed to be never a hypernucleus)
391 G4int anA=theNucleus->GetMassNumber();
392 G4int aZ=theNucleus->GetCharge();
393 //G4int aL=theNucleus->GetNumberOfLambdas(); // Should be 0
394 G4int numberOfEx = 0;
395 G4int numberOfCh = 0;
396 G4int numberOfHoles = 0;
397 G4double exEnergy = 0.0;
398 G4double R = theNucleus->GetNuclearRadius();
399 G4LorentzVector Target4Momentum(0.,0.,0.,0.);
400
401 // loop over the wounded target nucleus
402 G4Nucleon * theCurrentNucleon =
403 theNucleus->StartLoop() ? theNucleus->GetNextNucleon() : 0;
404 while(theCurrentNucleon) /* Loop checking, 31.08.2015, G.Folger */
405 {
406 if(theCurrentNucleon->AreYouHit()) {
407 ++numberOfHoles;
408 ++numberOfEx;
409 --anA;
410 aZ -= G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/
411 eplus + 0.1);
412 exEnergy += theCurrentNucleon->GetBindingEnergy();
413 Target4Momentum -=theCurrentNucleon->Get4Momentum();
414 }
415 theCurrentNucleon = theNucleus->GetNextNucleon();
416 }
417
418#ifdef debugPrecoInt
419 G4cout<<"Residual Target A Z (numberOfLambdas) E* 4mom "<<anA<<" "<<aZ<<" (0"//<<aL
420 <<") "<<exEnergy<<" "<<Target4Momentum<<G4endl;
421#endif
422
423 // prepare the projectile residual - which can be a hypernucleus or anti-hypernucleus
424
425 G4bool ProjectileIsAntiNucleus=
427
429
430 G4int anAb=theProjectileNucleus->GetMassNumber();
431 G4int aZb=theProjectileNucleus->GetCharge();
432 G4int aLb=theProjectileNucleus->GetNumberOfLambdas(); // Non negative number of (anti-)lambdas in (anti-)nucleus
433 G4int numberOfExB = 0;
434 G4int numberOfChB = 0;
435 G4int numberOfHolesB = 0;
436 G4double exEnergyB = 0.0;
437 G4double Rb = theProjectileNucleus->GetNuclearRadius();
438 G4LorentzVector Projectile4Momentum(0.,0.,0.,0.);
439
440 // loop over the wounded projectile nucleus or anti-nucleus
441 theCurrentNucleon =
442 theProjectileNucleus->StartLoop() ? theProjectileNucleus->GetNextNucleon() : 0;
443 while(theCurrentNucleon) /* Loop checking, 31.08.2015, G.Folger */
444 {
445 if(theCurrentNucleon->AreYouHit()) {
446 ++numberOfHolesB;
447 ++numberOfExB;
448 --anAb;
449 if(!ProjectileIsAntiNucleus) {
450 aZb -= G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/
451 eplus + 0.1);
452 if (theCurrentNucleon->GetParticleType()==G4Lambda::Definition()) --aLb;
453 } else {
454 aZb += G4int(theCurrentNucleon->GetDefinition()->GetPDGCharge()/
455 eplus - 0.1);
456 if (theCurrentNucleon->GetParticleType()==G4AntiLambda::Definition()) --aLb;
457 }
458 exEnergyB += theCurrentNucleon->GetBindingEnergy();
459 Projectile4Momentum -=theCurrentNucleon->Get4Momentum();
460 }
461 theCurrentNucleon = theProjectileNucleus->GetNextNucleon();
462 }
463
464 G4bool ExistTargetRemnant = G4double (numberOfHoles) <
465 0.3* G4double (numberOfHoles + anA);
466 G4bool ExistProjectileRemnant= G4double (numberOfHolesB) <
467 0.3*G4double (numberOfHolesB + anAb);
468
469#ifdef debugPrecoInt
470 G4cout<<"Projectile residual A Z (numberOfLambdas) E* 4mom "<<anAb<<" "<<aZb<<" ("<<aLb
471 <<") "<<exEnergyB<<" "<<Projectile4Momentum<<G4endl;
472 G4cout<<" ExistTargetRemnant ExistProjectileRemnant "
473 <<ExistTargetRemnant<<" "<< ExistProjectileRemnant<<G4endl;
474#endif
475 //-----------------------------------------------------------------------------
476 // decay the strong resonances
478 G4DecayKineticTracks decay(theSecondaries);
479
480 MakeCoalescence(theSecondaries);
481
482 #ifdef debugPrecoInt
483 G4cout<<"Secondary stable particles number "<<theSecondaries->size()<<G4endl;
484 #endif
485
486#ifdef debugPrecoInt
487 G4LorentzVector secondary4Momemtum(0,0,0,0);
488 G4int SecondrNum(0);
489#endif
490
491 // Loop over secondaries.
492 // We are assuming that only protons and neutrons - for nuclei -
493 // and only antiprotons and antineutrons - for antinuclei - can be absorbed,
494 // not instead lambdas (or hyperons more generally) - for nuclei - or anti-lambdas
495 // (or anti-hyperons more generally) - for antinuclei. This is a simplification,
496 // to be eventually reviewed later on, in particular when generic hypernuclei and
497 // anti-hypernuclei are introduced, instead of the few light hypernuclei and
498 // anti-hypernuclei which currently exist in Geant4.
499 G4KineticTrackVector::iterator iter;
500 for(iter=theSecondaries->begin(); iter !=theSecondaries->end(); ++iter)
501 {
502 const G4ParticleDefinition* part = (*iter)->GetDefinition();
503 G4LorentzVector aTrack4Momentum=(*iter)->Get4Momentum();
504
505 if( part != proton && part != neutron &&
506 (part != ANTIproton && ProjectileIsAntiNucleus) &&
507 (part != ANTIneutron && ProjectileIsAntiNucleus) )
508 {
509 G4ReactionProduct * theNew = new G4ReactionProduct(part);
510 theNew->SetMomentum(aTrack4Momentum.vect());
511 theNew->SetTotalEnergy(aTrack4Momentum.e());
512 theNew->SetCreatorModelID((*iter)->GetCreatorModelID());
513 theTotalResult->push_back(theNew);
514#ifdef debugPrecoInt
515 SecondrNum++;
516 secondary4Momemtum += (*iter)->Get4Momentum();
517 G4cout<<"Secondary "<<SecondrNum<<" "
518 <<theNew->GetDefinition()->GetParticleName()<<" "
519 <<theNew->GetMomentum()<<" "<<theNew->GetTotalEnergy()<<G4endl;
520
521#endif
522 delete (*iter);
523 continue;
524 }
525
526 G4bool CanBeCapturedByTarget = false;
527 if( part == proton || part == neutron)
528 {
529 CanBeCapturedByTarget = ExistTargetRemnant &&
531 (aTrack4Momentum + Target4Momentum).mag() -
532 aTrack4Momentum.mag() - Target4Momentum.mag()) &&
533 ((*iter)->GetPosition().mag() < R);
534 }
535 // ---------------------------
536 G4LorentzVector Position((*iter)->GetPosition(), (*iter)->GetFormationTime());
537 Position.boost(bst);
538
539 G4bool CanBeCapturedByProjectile = false;
540
541 if( !ProjectileIsAntiNucleus &&
542 ( part == proton || part == neutron))
543 {
544 CanBeCapturedByProjectile = ExistProjectileRemnant &&
546 (aTrack4Momentum + Projectile4Momentum).mag() -
547 aTrack4Momentum.mag() - Projectile4Momentum.mag()) &&
548 (Position.vect().mag() < Rb);
549 }
550
551 if( ProjectileIsAntiNucleus &&
552 ( part == ANTIproton || part == ANTIneutron))
553 {
554 CanBeCapturedByProjectile = ExistProjectileRemnant &&
556 (aTrack4Momentum + Projectile4Momentum).mag() -
557 aTrack4Momentum.mag() - Projectile4Momentum.mag()) &&
558 (Position.vect().mag() < Rb);
559 }
560
561 if(CanBeCapturedByTarget && CanBeCapturedByProjectile)
562 {
563 if(G4UniformRand() < 0.5)
564 { CanBeCapturedByTarget = true; CanBeCapturedByProjectile = false;}
565 else
566 { CanBeCapturedByTarget = false; CanBeCapturedByProjectile = true;}
567 }
568
569 if(CanBeCapturedByTarget)
570 {
571 // within the target nucleus, neutron or proton
572 // now calculate A, Z of the fragment, momentum,
573 // number of exciton states
574#ifdef debugPrecoInt
575 G4cout<<"Track is CapturedByTarget "<<" "<<part->GetParticleName()<<" "
576 <<aTrack4Momentum<<" "<<aTrack4Momentum.mag()<<G4endl;
577#endif
578 ++anA;
579 ++numberOfEx;
580 G4int Z = G4int(part->GetPDGCharge()/eplus + 0.1);
581 aZ += Z;
582 numberOfCh += Z;
583 Target4Momentum +=aTrack4Momentum;
584 delete (*iter);
585 } else if(CanBeCapturedByProjectile)
586 {
587 // within the projectile nucleus, neutron or proton
588 // now calculate A, Z of the fragment, momentum,
589 // number of exciton states
590#ifdef debugPrecoInt
591 G4cout<<"Track is CapturedByProjectile"<<" "<<part->GetParticleName()<<" "
592 <<aTrack4Momentum<<" "<<aTrack4Momentum.mag()<<G4endl;
593#endif
594 ++anAb;
595 ++numberOfExB;
596 G4int Z = G4int(part->GetPDGCharge()/eplus + 0.1);
597 if( ProjectileIsAntiNucleus ) Z=-Z;
598 aZb += Z;
599 numberOfChB += Z;
600 Projectile4Momentum +=aTrack4Momentum;
601 delete (*iter);
602 } else
603 { // the track is not captured
604 G4ReactionProduct * theNew = new G4ReactionProduct(part);
605 theNew->SetMomentum(aTrack4Momentum.vect());
606 theNew->SetTotalEnergy(aTrack4Momentum.e());
607 theNew->SetCreatorModelID((*iter)->GetCreatorModelID());
608 theTotalResult->push_back(theNew);
609
610#ifdef debugPrecoInt
611 SecondrNum++;
612 secondary4Momemtum += (*iter)->Get4Momentum();
613/*
614 G4cout<<"Secondary "<<SecondrNum<<" "
615 <<theNew->GetDefinition()->GetParticleName()<<" "
616 <<secondary4Momemtum<<G4endl;
617*/
618#endif
619 delete (*iter);
620 continue;
621 }
622 }
623 delete theSecondaries;
624 //-----------------------------------------------------
625
626 #ifdef debugPrecoInt
627 G4cout<<"Final target residual A Z (numberOfLambdas) E* 4mom "<<anA<<" "<<aZ<<" (0"//<<aL
628 <<") "<<exEnergy<<" "<<Target4Momentum<<G4endl;
629 #endif
630
631 if(0!=anA )
632 {
633 // We assume that the target residual is never a hypernucleus
635
636 if((anA == theNucleus->GetMassNumber()) && (exEnergy <= 0.))
637 {Target4Momentum.setE(fMass);}
638
639 G4double RemnMass=Target4Momentum.mag();
640
641 if(RemnMass < fMass)
642 {
643 RemnMass=fMass + exEnergy;
644 Target4Momentum.setE(std::sqrt(Target4Momentum.vect().mag2() +
645 RemnMass*RemnMass));
646 } else
647 { exEnergy=RemnMass-fMass;}
648
649 if( exEnergy < 0.) exEnergy=0.;
650
651 // Need to de-excite the remnant nucleus
652 G4Fragment anInitialState(anA, aZ, Target4Momentum);
653 anInitialState.SetNumberOfParticles(numberOfEx-numberOfHoles);
654 anInitialState.SetNumberOfCharged(numberOfCh);
655 anInitialState.SetNumberOfHoles(numberOfHoles);
656 anInitialState.SetCreatorModelID(secID);
657
658 G4ReactionProductVector * aPrecoResult =
659 theDeExcitation->DeExcite(anInitialState);
660
661 #ifdef debugPrecoInt
662 G4cout<<"Target fragment number "<<aPrecoResult->size()<<G4endl;
663 #endif
664
665 // fill pre-compound part into the result, and return
666 for(unsigned int ll=0; ll<aPrecoResult->size(); ++ll)
667 {
668 theTotalResult->push_back(aPrecoResult->operator[](ll));
669 #ifdef debugPrecoInt
670 G4cout<<"Target fragment "<<ll<<" "
671 <<aPrecoResult->operator[](ll)->GetDefinition()->GetParticleName()<<" "
672 <<aPrecoResult->operator[](ll)->GetMomentum()<<" "
673 <<aPrecoResult->operator[](ll)->GetTotalEnergy()<<" "
674 <<aPrecoResult->operator[](ll)->GetMass()<<G4endl;
675 #endif
676 }
677 delete aPrecoResult;
678 }
679
680 //-----------------------------------------------------
681 if((anAb == theProjectileNucleus->GetMassNumber())&& (exEnergyB <= 0.))
682 {Projectile4Momentum = GetPrimaryProjectile()->Get4Momentum();}
683
684 #ifdef debugPrecoInt
685 G4cout<<"Final projectile residual A Z (numberOfLambdas) E* Pmom Pmag2 "<<anAb<<" "<<aZb<<" ("
686 <<aLb<<") "<<exEnergyB<<" "<<Projectile4Momentum<<" "
687 <<Projectile4Momentum.mag2()<<G4endl;
688 #endif
689
690 if(0!=anAb)
691 {
692 // The projectile residual can be a hypernucleus or anti-hypernucleus
693 G4double fMass = 0.0;
694 if ( aLb > 0 ) {
695 fMass = G4HyperNucleiProperties::GetNuclearMass(anAb, aZb, aLb);
696 } else {
697 fMass = G4NucleiProperties::GetNuclearMass(anAb, aZb);
698 }
699 G4double RemnMass=Projectile4Momentum.mag();
700
701 if(RemnMass < fMass)
702 {
703 RemnMass=fMass + exEnergyB;
704 Projectile4Momentum.setE(std::sqrt(Projectile4Momentum.vect().mag2() +
705 RemnMass*RemnMass));
706 } else
707 { exEnergyB=RemnMass-fMass;}
708
709 if( exEnergyB < 0.) exEnergyB=0.;
710
711 G4ThreeVector bstToCM =Projectile4Momentum.findBoostToCM();
712 Projectile4Momentum.boost(bstToCM);
713
714 // Need to de-excite the remnant nucleus
715 G4Fragment anInitialState(anAb, aZb, aLb, Projectile4Momentum);
716 anInitialState.SetNumberOfParticles(numberOfExB-numberOfHolesB);
717 anInitialState.SetNumberOfCharged(numberOfChB);
718 anInitialState.SetNumberOfHoles(numberOfHolesB);
719 anInitialState.SetCreatorModelID(secID);
720
721 G4ReactionProductVector * aPrecoResult =
722 theDeExcitation->DeExcite(anInitialState);
723
724 #ifdef debugPrecoInt
725 G4cout<<"Projectile fragment number "<<aPrecoResult->size()<<G4endl;
726 #endif
727
728 // fill pre-compound part into the result, and return
729 for(unsigned int ll=0; ll<aPrecoResult->size(); ++ll)
730 {
731 G4LorentzVector tmp=G4LorentzVector(aPrecoResult->operator[](ll)->GetMomentum(),
732 aPrecoResult->operator[](ll)->GetTotalEnergy());
733 tmp.boost(-bstToCM); // Transformation to the system of original remnant
734 aPrecoResult->operator[](ll)->SetMomentum(tmp.vect());
735 aPrecoResult->operator[](ll)->SetTotalEnergy(tmp.e());
736
737 if(ProjectileIsAntiNucleus)
738 {
739 const G4ParticleDefinition * aFragment=aPrecoResult->operator[](ll)->GetDefinition();
740 const G4ParticleDefinition * LastFragment=aFragment;
741 if (aFragment == proton) {LastFragment=G4AntiProton::AntiProtonDefinition();}
742 else if(aFragment == neutron) {LastFragment=G4AntiNeutron::AntiNeutronDefinition();}
743 else if(aFragment == deuteron){LastFragment=G4AntiDeuteron::AntiDeuteronDefinition();}
744 else if(aFragment == triton) {LastFragment=G4AntiTriton::AntiTritonDefinition();}
745 else if(aFragment == He3) {LastFragment=G4AntiHe3::AntiHe3Definition();}
746 else if(aFragment == He4) {LastFragment=G4AntiAlpha::AntiAlphaDefinition();}
747 else {}
748
749 if (aLb > 0) { // Anti-hypernucleus
750 if (aFragment == G4HyperTriton::Definition()) {
751 LastFragment=G4AntiHyperTriton::Definition();
752 } else if (aFragment == G4HyperH4::Definition()) {
753 LastFragment=G4AntiHyperH4::Definition();
754 } else if (aFragment == G4HyperAlpha::Definition()) {
755 LastFragment=G4AntiHyperAlpha::Definition();
756 } else if (aFragment == G4HyperHe5::Definition()) {
757 LastFragment=G4AntiHyperHe5::Definition();
758 } else if (aFragment == G4DoubleHyperH4::Definition()) {
759 LastFragment=G4AntiDoubleHyperH4::Definition();
760 } else if (aFragment == G4DoubleHyperDoubleNeutron::Definition()) {
762 }
763 }
764
765 aPrecoResult->operator[](ll)->SetDefinitionAndUpdateE(LastFragment);
766 }
767
768 #ifdef debugPrecoInt
769 G4cout<<"Projectile fragment "<<ll<<" "
770 <<aPrecoResult->operator[](ll)->GetDefinition()->GetParticleName()<<" "
771 <<aPrecoResult->operator[](ll)->GetMomentum()<<" "
772 <<aPrecoResult->operator[](ll)->GetTotalEnergy()<<" "
773 <<aPrecoResult->operator[](ll)->GetMass()<<G4endl;
774 #endif
775
776 theTotalResult->push_back(aPrecoResult->operator[](ll));
777 }
778
779 delete aPrecoResult;
780 }
781
782 return theTotalResult;
783}
784
785
787 // This method replaces pairs of proton-neutron - in the case of nuclei - or
788 // antiproton-antineutron - in the case of anti-nuclei - which are close in
789 // momentum, with, respectively, deuterons and anti-deuterons.
790 // Note that in the case of hypernuclei or anti-hypernuclei, lambdas or anti-lambdas
791 // are not considered for coalescence because hyper-deuteron or anti-hyper-deuteron
792 // are assumed not to exist.
793
794 if (!tracks) return;
795
796 G4double MassCut = deuteron->GetPDGMass() + DeltaM; // In MeV
797
798 for ( size_t i = 0; i < tracks->size(); ++i ) { // search for protons
799
800 G4KineticTrack* trackP = (*tracks)[i];
801 if ( ! trackP ) continue;
802 if (trackP->GetDefinition() != proton) continue;
803
804 G4LorentzVector Prot4Mom = trackP->Get4Momentum();
805 G4LorentzVector ProtSPposition = G4LorentzVector(trackP->GetPosition(), trackP->GetFormationTime());
806
807 for ( size_t j = 0; j < tracks->size(); ++j ) { // search for neutron
808
809 G4KineticTrack* trackN = (*tracks)[j];
810 if (! trackN ) continue;
811 if (trackN->GetDefinition() != neutron) continue;
812
813 G4LorentzVector Neut4Mom = trackN->Get4Momentum();
814 G4LorentzVector NeutSPposition = G4LorentzVector( trackN->GetPosition(), trackN->GetFormationTime()*hbarc/fermi);
815 G4double EffMass = (Prot4Mom + Neut4Mom).mag();
816
817 if ( EffMass <= MassCut ) { // && (EffDistance <= SpaceCut)) { // Create deuteron
818 G4KineticTrack* aDeuteron =
820 (trackP->GetFormationTime() + trackN->GetFormationTime())/2.0,
821 (trackP->GetPosition() + trackN->GetPosition() )/2.0,
822 ( Prot4Mom + Neut4Mom ));
823 aDeuteron->SetCreatorModelID(secID);
824 tracks->push_back(aDeuteron);
825 delete trackP; delete trackN;
826 (*tracks)[i] = nullptr; (*tracks)[j] = nullptr;
827 break;
828 }
829 }
830 }
831
832 // Find and remove null pointers created by decays above
833 for ( int jj = tracks->size()-1; jj >= 0; --jj ) {
834 if ( ! (*tracks)[jj] ) tracks->erase(tracks->begin()+jj);
835 }
836
837}
838
G4double G4Log(G4double x)
Definition: G4Log.hh:226
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double fermi
Definition: G4SIunits.hh:83
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4cin
Definition: G4ios.hh:56
#define G4UniformRand()
Definition: Randomize.hh:52
double mag2() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
Hep3Vector findBoostToCM() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4AntiAlpha * AntiAlpha()
Definition: G4AntiAlpha.cc:88
static G4AntiAlpha * AntiAlphaDefinition()
Definition: G4AntiAlpha.cc:83
static G4AntiDeuteron * AntiDeuteronDefinition()
static G4AntiDeuteron * AntiDeuteron()
static G4AntiDoubleHyperDoubleNeutron * Definition()
static G4AntiDoubleHyperH4 * Definition()
static G4AntiHe3 * AntiHe3()
Definition: G4AntiHe3.cc:93
static G4AntiHe3 * AntiHe3Definition()
Definition: G4AntiHe3.cc:88
static G4AntiHyperAlpha * Definition()
static G4AntiHyperH4 * Definition()
static G4AntiHyperHe5 * Definition()
static G4AntiHyperTriton * Definition()
static G4AntiLambda * Definition()
Definition: G4AntiLambda.cc:52
static G4AntiNeutron * AntiNeutronDefinition()
static G4AntiNeutron * AntiNeutron()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:92
static G4AntiProton * AntiProtonDefinition()
Definition: G4AntiProton.cc:87
static G4AntiTriton * AntiTriton()
Definition: G4AntiTriton.cc:93
static G4AntiTriton * AntiTritonDefinition()
Definition: G4AntiTriton.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4DoubleHyperDoubleNeutron * Definition()
static G4DoubleHyperH4 * Definition()
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:405
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:391
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:400
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus)
virtual G4ReactionProductVector * PropagateNuclNucl(G4KineticTrackVector *theSecondaries, G4V3DNucleus *theNucleus, G4V3DNucleus *theProjectileNucleus)
virtual void PropagateModelDescription(std::ostream &) const
G4GeneratorPrecompoundInterface(G4VPreCompoundModel *p=0)
void MakeCoalescence(G4KineticTrackVector *theSecondaries)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
static G4HyperAlpha * Definition()
Definition: G4HyperAlpha.cc:45
static G4HyperH4 * Definition()
Definition: G4HyperH4.cc:45
static G4HyperHe5 * Definition()
Definition: G4HyperHe5.cc:45
static G4double GetNuclearMass(G4int A, G4int Z, G4int L)
static G4HyperTriton * Definition()
G4double GetFormationTime() const
void SetCreatorModelID(G4int id)
const G4ThreeVector & GetPosition() const
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
static G4Lambda * Definition()
Definition: G4Lambda.cc:52
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4bool AreYouHit() const
Definition: G4Nucleon.hh:98
const G4ParticleDefinition * GetParticleType() const
Definition: G4Nucleon.hh:85
virtual const G4LorentzVector & Get4Momentum() const
Definition: G4Nucleon.hh:72
G4double GetBindingEnergy() const
Definition: G4Nucleon.hh:75
virtual const G4ParticleDefinition * GetDefinition() const
Definition: G4Nucleon.hh:86
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
void SetCreatorModelID(const G4int mod)
G4ThreeVector GetMomentum() const
static G4Triton * Triton()
Definition: G4Triton.cc:93
virtual G4Nucleon * GetNextNucleon()=0
virtual G4int GetCharge()=0
virtual G4int GetNumberOfLambdas()=0
virtual G4bool StartLoop()=0
virtual G4double GetNuclearRadius()=0
virtual G4int GetMassNumber()=0
const G4HadProjectile * GetPrimaryProjectile() const
void SetDeExcitation(G4VPreCompoundModel *ptr)
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)=0
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
float hbarc
Definition: hepunit.py:264
T sqr(const T &x)
Definition: templates.hh:128