Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QMDGroundStateNucleus.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// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
27//
29
30#include "G4NucleiProperties.hh"
31#include "G4Proton.hh"
32#include "G4Neutron.hh"
33#include "G4Exp.hh"
34#include "G4Pow.hh"
36#include "Randomize.hh"
37
39: maxTrial ( 1000 )
40, r00 ( 1.124 ) // radius parameter for Woods-Saxon [fm]
41, r01 ( 0.5 ) // radius parameter for Woods-Saxon
42, saa ( 0.2 ) // diffuse parameter for initial Woods-Saxon shape
43, rada ( 0.9 ) // cutoff parameter
44, radb ( 0.3 ) // cutoff parameter
45, dsam ( 1.5 ) // minimum distance for same particle [fm]
46, ddif ( 1.0 ) // minimum distance for different particle
47, edepth ( 0.0 )
48, epse ( 0.000001 ) // torelance for energy in [GeV]
49, meanfield ( NULL )
50{
51
52 //std::cout << " G4QMDGroundStateNucleus( G4int z , G4int a ) Begin " << z << " " << a << std::endl;
53
54 dsam2 = dsam*dsam;
55 ddif2 = ddif*ddif;
56
58
59 hbc = parameters->Get_hbc();
60 gamm = parameters->Get_gamm();
61 cpw = parameters->Get_cpw();
62 cph = parameters->Get_cph();
63 epsx = parameters->Get_epsx();
64 cpc = parameters->Get_cpc();
65
66 cdp = parameters->Get_cdp();
67 c0p = parameters->Get_c0p();
68 c3p = parameters->Get_c3p();
69 csp = parameters->Get_csp();
70 clp = parameters->Get_clp();
71
72 // Following 10 lines should be here, right before the line 90.
73 // Otherwise, mass number cannot be conserved if the projectile or
74 // the target are nucleons.
75 //Nucleon primary or target case;
76 if ( z == 1 && a == 1 ) { // Hydrogen Case or proton primary
78 ebini = 0.0;
79 return;
80 }
81 else if ( z == 0 && a == 1 ) { // Neutron primary
83 ebini = 0.0;
84 return;
85 }
86
87
88 //edepth = 0.0;
89
90 for ( int i = 0 ; i < a ; i++ )
91 {
92
94
95 if ( i < z )
96 {
97 pd = G4Proton::Proton();
98 }
99 else
100 {
101 pd = G4Neutron::Neutron();
102 }
103
104 G4ThreeVector p( 0.0 );
105 G4ThreeVector r( 0.0 );
106 G4QMDParticipant* aParticipant = new G4QMDParticipant( pd , p , r );
107 SetParticipant( aParticipant );
108
109 }
110
111 G4double radious = r00 * G4Pow::GetInstance()->A13( double ( GetMassNumber() ) );
112
113 rt00 = radious - r01;
114 radm = radious - rada * ( gamm - 1.0 ) + radb;
115 rmax = 1.0 / ( 1.0 + G4Exp ( -rt00/saa ) );
116
117 //maxTrial = 1000;
118
119
120 meanfield = new G4QMDMeanField();
121 meanfield->SetSystem( this );
122
123 //std::cout << "G4QMDGroundStateNucleus( G4int z , G4int a ) packNucleons Begin ( z , a ) ( " << z << ", " << a << " )" << std::endl;
124 packNucleons();
125 //std::cout << "G4QMDGroundStateNucleus( G4int z , G4int a ) packNucleons End" << std::endl;
126
127 delete meanfield;
128
129}
130
131
132
133void G4QMDGroundStateNucleus::packNucleons()
134{
135
136 //std::cout << "G4QMDGroundStateNucleus::packNucleons" << std::endl;
137
139
140 G4double ebin00 = ebini * 0.001;
141
142 G4double ebin0 = 0.0;
143 G4double ebin1 = 0.0;
144
145 if ( GetMassNumber() != 4 )
146 {
147 ebin0 = ( ebini - 0.5 ) * 0.001;
148 ebin1 = ( ebini + 0.5 ) * 0.001;
149 }
150 else
151 {
152 ebin0 = ( ebini - 1.5 ) * 0.001;
153 ebin1 = ( ebini + 1.5 ) * 0.001;
154 }
155
156 G4int n0Try = 0;
157 G4bool isThisOK = false;
158 while ( n0Try < maxTrial ) // Loop checking, 11.03.2015, T. Koi
159 {
160 n0Try++;
161 //std::cout << "TKDB packNucleons n0Try " << n0Try << std::endl;
162
163// Sampling Position
164
165 //std::cout << "TKDB Sampling Position " << std::endl;
166
167 G4bool areThesePsOK = false;
168 G4int npTry = 0;
169 while ( npTry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
170 {
171 //std::cout << "TKDB Sampling Position npTry " << npTry << std::endl;
172 npTry++;
173 G4int i = 0;
174 if ( samplingPosition( i ) )
175 {
176 //std::cout << "packNucleons samplingPosition 0 succeed " << std::endl;
177 for ( i = 1 ; i < GetMassNumber() ; i++ )
178 {
179 //std::cout << "packNucleons samplingPosition " << i << " trying " << std::endl;
180 if ( !( samplingPosition( i ) ) )
181 {
182 //std::cout << "packNucleons samplingPosition " << i << " failed" << std::endl;
183 break;
184 }
185 }
186 if ( i == GetMassNumber() )
187 {
188 //std::cout << "packNucleons samplingPosition all scucceed " << std::endl;
189 areThesePsOK = true;
190 break;
191 }
192 }
193 }
194 if ( areThesePsOK == false ) continue;
195
196 //std::cout << "TKDB Sampling Position End" << std::endl;
197
198// Calculate Two-body quantities
199
200 meanfield->Cal2BodyQuantities();
201 std::vector< G4double > rho_a ( GetMassNumber() , 0.0 );
202 std::vector< G4double > rho_s ( GetMassNumber() , 0.0 );
203 std::vector< G4double > rho_c ( GetMassNumber() , 0.0 );
204
205 rho_l.resize ( GetMassNumber() , 0.0 );
206 d_pot.resize ( GetMassNumber() , 0.0 );
207
208 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
209 {
210 for ( G4int j = 0 ; j < GetMassNumber() ; j++ )
211 {
212
213 rho_a[ i ] += meanfield->GetRHA( j , i );
214 G4int k = 0;
215 if ( participants[j]->GetDefinition() != participants[i]->GetDefinition() )
216 {
217 k = 1;
218 }
219
220 rho_s[ i ] += meanfield->GetRHA( j , i )*( 1.0 - 2.0 * k ); // OK?
221
222 rho_c[ i ] += meanfield->GetRHE( j , i );
223 }
224
225 }
226
227 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
228 {
229 rho_l[i] = cdp * ( rho_a[i] + 1.0 );
230 d_pot[i] = c0p * rho_a[i]
231 + c3p * G4Pow::GetInstance()->powA ( rho_a[i] , gamm )
232 + csp * rho_s[i]
233 + clp * rho_c[i];
234
235 //std::cout << "d_pot[i] " << i << " " << d_pot[i] << std::endl;
236 }
237
238
239// Sampling Momentum
240
241 //std::cout << "TKDB Sampling Momentum " << std::endl;
242
243 phase_g.clear();
244 phase_g.resize( GetMassNumber() , 0.0 );
245
246 //std::cout << "TKDB Sampling Momentum 1st " << std::endl;
247 G4bool isThis1stMOK = false;
248 G4int nmTry = 0;
249 while ( nmTry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
250 {
251 nmTry++;
252 G4int i = 0;
253 if ( samplingMomentum( i ) )
254 {
255 isThis1stMOK = true;
256 break;
257 }
258 }
259 if ( isThis1stMOK == false ) continue;
260
261 //std::cout << "TKDB Sampling Momentum 2nd so on" << std::endl;
262
263 G4bool areTheseMsOK = true;
264 nmTry = 0;
265 while ( nmTry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
266 {
267 nmTry++;
268 G4int i = 0;
269 for ( i = 1 ; i < GetMassNumber() ; i++ )
270 {
271 //std::cout << "TKDB packNucleons samplingMomentum try " << i << std::endl;
272 if ( !( samplingMomentum( i ) ) )
273 {
274 //std::cout << "TKDB packNucleons samplingMomentum " << i << " failed" << std::endl;
275 areTheseMsOK = false;
276 break;
277 }
278 }
279 if ( i == GetMassNumber() )
280 {
281 areTheseMsOK = true;
282 }
283
284 break;
285 }
286 if ( areTheseMsOK == false ) continue;
287
288// Kill Angluar Momentum
289
290 //std::cout << "TKDB Sampling Kill Angluar Momentum " << std::endl;
291
292 killCMMotionAndAngularM();
293
294
295// Check Binding Energy
296
297 //std::cout << "packNucleons Check Binding Energy Begin " << std::endl;
298
299 G4double ekinal = 0.0;
300 for ( int i = 0 ; i < GetMassNumber() ; i++ )
301 {
302 ekinal += participants[i]->GetKineticEnergy();
303 }
304
305 meanfield->Cal2BodyQuantities();
306
307 G4double totalPotentialE = meanfield->GetTotalPotential();
308 G4double ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
309
310 //std::cout << "packNucleons totalPotentialE " << totalPotentialE << std::endl;
311 //std::cout << "packNucleons ebinal " << ebinal << std::endl;
312 //std::cout << "packNucleons ekinal " << ekinal << std::endl;
313
314 if ( ebinal < ebin0 || ebinal > ebin1 )
315 {
316 //std::cout << "packNucleons ebin0 " << ebin0 << std::endl;
317 //std::cout << "packNucleons ebin1 " << ebin1 << std::endl;
318 //std::cout << "packNucleons ebinal " << ebinal << std::endl;
319 //std::cout << "packNucleons Check Binding Energy Failed " << std::endl;
320 continue;
321 }
322
323 //std::cout << "packNucleons Check Binding Energy End = OK " << std::endl;
324
325
326// Energy Adujstment
327
328 G4double dtc = 1.0;
329 G4double frg = -0.1;
330 G4double rdf0 = 0.5;
331
332 G4double edif0 = ebinal - ebin00;
333
334 G4double cfrc = 0.0;
335 if ( 0 < edif0 )
336 {
337 cfrc = frg;
338 }
339 else
340 {
341 cfrc = -frg;
342 }
343
344 G4int ifrc = 1;
345
346 G4int neaTry = 0;
347
348 G4bool isThisEAOK = false;
349 while ( neaTry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
350 {
351 neaTry++;
352
353 G4double edif = ebinal - ebin00;
354
355 //std::cout << "TKDB edif " << edif << std::endl;
356 if ( std::abs ( edif ) < epse )
357 {
358
359 isThisEAOK = true;
360 //std::cout << "isThisEAOK " << isThisEAOK << std::endl;
361 break;
362 }
363
364 G4int jfrc = 0;
365 if ( edif < 0.0 )
366 {
367 jfrc = 1;
368 }
369 else
370 {
371 jfrc = -1;
372 }
373
374 if ( jfrc != ifrc )
375 {
376 cfrc = -rdf0 * cfrc;
377 dtc = rdf0 * dtc;
378 }
379
380 if ( jfrc == ifrc && std::abs( edif0 ) < std::abs( edif ) )
381 {
382 cfrc = -rdf0 * cfrc;
383 dtc = rdf0 * dtc;
384 }
385
386 ifrc = jfrc;
387 edif0 = edif;
388
389 meanfield->CalGraduate();
390
391 for ( int i = 0 ; i < GetMassNumber() ; i++ )
392 {
393 G4ThreeVector ri = participants[i]->GetPosition();
394 G4ThreeVector p_i = participants[i]->GetMomentum();
395
396 ri += dtc * ( meanfield->GetFFr(i) - cfrc * ( meanfield->GetFFp(i) ) );
397 p_i += dtc * ( meanfield->GetFFp(i) + cfrc * ( meanfield->GetFFr(i) ) );
398
399 participants[i]->SetPosition( ri );
400 participants[i]->SetMomentum( p_i );
401 }
402
403 ekinal = 0.0;
404
405 for ( int i = 0 ; i < GetMassNumber() ; i++ )
406 {
407 ekinal += participants[i]->GetKineticEnergy();
408 }
409
410 meanfield->Cal2BodyQuantities();
411 totalPotentialE = meanfield->GetTotalPotential();
412
413 ebinal = ( totalPotentialE + ekinal ) / double ( GetMassNumber() );
414
415 }
416 //std::cout << "isThisEAOK " << isThisEAOK << std::endl;
417 if ( isThisEAOK == false ) continue;
418
419 isThisOK = true;
420 //std::cout << "isThisOK " << isThisOK << std::endl;
421 break;
422
423 }
424
425 if ( isThisOK == false )
426 {
427 G4cout << "GroundStateNucleus state cannot be created. Try again with another parameters." << G4endl;
428 }
429
430 //std::cout << "packNucleons End" << std::endl;
431 return;
432}
433
434
435G4bool G4QMDGroundStateNucleus::samplingPosition( G4int i )
436{
437
438 G4bool result = false;
439
440 G4int nTry = 0;
441
442 while ( nTry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
443 {
444
445 //std::cout << "samplingPosition i th particle, nTtry " << i << " " << nTry << std::endl;
446
447 G4double rwod = -1.0;
448 G4double rrr = 0.0;
449
450 G4double rx = 0.0;
451 G4double ry = 0.0;
452 G4double rz = 0.0;
453
454 G4int icounter = 0;
455 G4int icounter_max = 1024;
456 while ( G4UniformRand() * rmax > rwod ) // Loop checking, 11.03.2015, T. Koi
457 {
458 icounter++;
459 if ( icounter > icounter_max ) {
460 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
461 break;
462 }
463
464 G4double rsqr = 10.0;
465 G4int jcounter = 0;
466 G4int jcounter_max = 1024;
467 while ( rsqr > 1.0 ) // Loop checking, 11.03.2015, T. Koi
468 {
469 jcounter++;
470 if ( jcounter > jcounter_max ) {
471 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
472 break;
473 }
474 rx = 1.0 - 2.0 * G4UniformRand();
475 ry = 1.0 - 2.0 * G4UniformRand();
476 rz = 1.0 - 2.0 * G4UniformRand();
477 rsqr = rx*rx + ry*ry + rz*rz;
478 }
479 rrr = radm * std::sqrt ( rsqr );
480 rwod = 1.0 / ( 1.0 + G4Exp ( ( rrr - rt00 ) / saa ) );
481
482 }
483
484 participants[i]->SetPosition( G4ThreeVector( rx , ry , rz )*radm );
485
486 if ( i == 0 )
487 {
488 result = true;
489 return result;
490 }
491
492// i > 1 ( Second Particle or later )
493// Check Distance to others
494
495 G4bool isThisOK = true;
496 for ( G4int j = 0 ; j < i ; j++ )
497 {
498
499 G4double r2 = participants[j]->GetPosition().diff2( participants[i]->GetPosition() );
500 G4double dmin2 = 0.0;
501
502 if ( participants[j]->GetDefinition() == participants[i]->GetDefinition() )
503 {
504 dmin2 = dsam2;
505 }
506 else
507 {
508 dmin2 = ddif2;
509 }
510
511 //std::cout << "distance between j and i " << j << " " << i << " " << r2 << " " << dmin2 << std::endl;
512 if ( r2 < dmin2 )
513 {
514 isThisOK = false;
515 break;
516 }
517
518 }
519
520 if ( isThisOK == true )
521 {
522 result = true;
523 return result;
524 }
525
526 nTry++;
527
528 }
529
530// Here return "false"
531 return result;
532}
533
534
535
536G4bool G4QMDGroundStateNucleus::samplingMomentum( G4int i )
537{
538
539 //std::cout << "TKDB samplingMomentum for " << i << std::endl;
540
541 G4bool result = false;
542
543 G4double pfm = hbc * G4Pow::GetInstance()->A13 ( ( 3.0 / 2.0 * pi*pi * rho_l[i] ) );
544
545 if ( 10 < GetMassNumber() && -5.5 < ebini )
546 {
547 pfm = pfm * ( 1.0 + 0.2 * std::sqrt( std::abs( 8.0 + ebini ) / 8.0 ) );
548 }
549
550 //std::cout << "TKDB samplingMomentum pfm " << pfm << std::endl;
551
552 std::vector< G4double > phase;
553 phase.resize( i+1 ); // i start from 0
554
555 G4int ntry = 0;
556// 710
557 while ( ntry < maxTrial ) // Loop checking, 11.03.2015, T. Koi
558 {
559
560 //std::cout << " TKDB ntry " << ntry << std::endl;
561 ntry++;
562
563 G4double ke = DBL_MAX;
564
565 G4int tkdb_i =0;
566// 700
567 G4int icounter = 0;
568 G4int icounter_max = 1024;
569 while ( ke + d_pot [i] > edepth ) // Loop checking, 11.03.2015, T. Koi
570 {
571 icounter++;
572 if ( icounter > icounter_max ) {
573 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
574 break;
575 }
576
577 G4double psqr = 10.0;
578 G4double px = 0.0;
579 G4double py = 0.0;
580 G4double pz = 0.0;
581
582 G4int jcounter = 0;
583 G4int jcounter_max = 1024;
584 while ( psqr > 1.0 ) // Loop checking, 11.03.2015, T. Koi
585 {
586 jcounter++;
587 if ( jcounter > jcounter_max ) {
588 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
589 break;
590 }
591 px = 1.0 - 2.0*G4UniformRand();
592 py = 1.0 - 2.0*G4UniformRand();
593 pz = 1.0 - 2.0*G4UniformRand();
594
595 psqr = px*px + py*py + pz*pz;
596 }
597
598 G4ThreeVector p ( px , py , pz );
599 p = pfm * p;
600 participants[i]->SetMomentum( p );
601 G4LorentzVector p4 = participants[i]->Get4Momentum();
602 //ke = p4.e() - p4.restMass();
603 ke = participants[i]->GetKineticEnergy();
604
605
606 tkdb_i++;
607 if ( tkdb_i > maxTrial ) return result; // return false
608
609 }
610
611 //std::cout << "TKDB ke d_pot[i] " << ke << " " << d_pot[i] << std::endl;
612
613
614 if ( i == 0 )
615 {
616 result = true;
617 return result;
618 }
619
620 G4bool isThisOK = true;
621
622 // Check Pauli principle
623
624 phase[ i ] = 0.0;
625
626 //std::cout << "TKDB Check Pauli Principle " << i << std::endl;
627
628 for ( G4int j = 0 ; j < i ; j++ )
629 {
630 phase[ j ] = 0.0;
631 //std::cout << "TKDB Check Pauli Principle i , j " << i << " , " << j << std::endl;
632 G4double expa = 0.0;
633 if ( participants[j]->GetDefinition() == participants[i]->GetDefinition() )
634 {
635
636 expa = - meanfield->GetRR2(i,j) * cpw;
637
638 if ( expa > epsx )
639 {
640 G4ThreeVector p_i = participants[i]->GetMomentum();
641 G4ThreeVector pj = participants[j]->GetMomentum();
642 G4double dist2_p = p_i.diff2( pj );
643
644 dist2_p = dist2_p*cph;
645 expa = expa - dist2_p;
646
647 if ( expa > epsx )
648 {
649
650 phase[j] = G4Exp ( expa );
651
652 if ( phase[j] * cpc > 0.2 )
653 {
654/*
655 G4cout << "TKDB Check Pauli Principle A i , j " << i << " , " << j << G4endl;
656 G4cout << "TKDB Check Pauli Principle phase[j] " << phase[j] << G4endl;
657 G4cout << "TKDB Check Pauli Principle phase[j]*cpc > 0.2 " << phase[j]*cpc << G4endl;
658*/
659 isThisOK = false;
660 break;
661 }
662 if ( ( phase_g[j] + phase[j] ) * cpc > 0.5 )
663 {
664/*
665 G4cout << "TKDB Check Pauli Principle B i , j " << i << " , " << j << G4endl;
666 G4cout << "TKDB Check Pauli Principle B phase_g[j] " << phase_g[j] << G4endl;
667 G4cout << "TKDB Check Pauli Principle B phase[j] " << phase[j] << G4endl;
668 G4cout << "TKDB Check Pauli Principle B phase_g[j] + phase[j] ) * cpc > 0.5 " << ( phase_g[j] + phase[j] ) * cpc << G4endl;
669*/
670 isThisOK = false;
671 break;
672 }
673
674 phase[i] += phase[j];
675 if ( phase[i] * cpc > 0.3 )
676 {
677/*
678 G4cout << "TKDB Check Pauli Principle C i , j " << i << " , " << j << G4endl;
679 G4cout << "TKDB Check Pauli Principle C phase[i] " << phase[i] << G4endl;
680 G4cout << "TKDB Check Pauli Principle C phase[i] * cpc > 0.3 " << phase[i] * cpc << G4endl;
681*/
682 isThisOK = false;
683 break;
684 }
685
686 //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
687
688 }
689 else
690 {
691 //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
692 }
693
694 }
695 else
696 {
697 //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
698 }
699
700 }
701 else
702 {
703 //std::cout << "TKDB Check Pauli Principle OK i , j " << i << " , " << j << std::endl;
704 }
705
706 }
707
708 if ( isThisOK == true )
709 {
710
711 phase_g[i] = phase[i];
712
713 for ( int j = 0 ; j < i ; j++ )
714 {
715 phase_g[j] += phase[j];
716 }
717
718 result = true;
719 return result;
720 }
721
722 }
723
724 return result;
725
726}
727
728
729
730void G4QMDGroundStateNucleus::killCMMotionAndAngularM()
731{
732
733// CalEnergyAndAngularMomentumInCM();
734
735 //std::vector< G4ThreeVector > p ( GetMassNumber() , 0.0 );
736 //std::vector< G4ThreeVector > r ( GetMassNumber() , 0.0 );
737
738// Move to cm system
739
740 G4ThreeVector pcm_tmp ( 0.0 );
741 G4ThreeVector rcm_tmp ( 0.0 );
742 G4double sumMass = 0.0;
743
744 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
745 {
746 pcm_tmp += participants[i]->GetMomentum();
747 rcm_tmp += participants[i]->GetPosition() * participants[i]->GetMass();
748 sumMass += participants[i]->GetMass();
749 }
750
751 pcm_tmp = pcm_tmp/GetMassNumber();
752 rcm_tmp = rcm_tmp/sumMass;
753
754 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
755 {
756 participants[i]->SetMomentum( participants[i]->GetMomentum() - pcm_tmp );
757 participants[i]->SetPosition( participants[i]->GetPosition() - rcm_tmp );
758 }
759
760// kill the angular momentum
761
762 G4ThreeVector ll ( 0.0 );
763 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
764 {
765 ll += participants[i]->GetPosition().cross ( participants[i]->GetMomentum() );
766 }
767
768 G4double rr[3][3];
769 G4double ss[3][3];
770 for ( G4int i = 0 ; i < 3 ; i++ )
771 {
772 for ( G4int j = 0 ; j < 3 ; j++ )
773 {
774 rr [i][j] = 0.0;
775
776 if ( i == j )
777 {
778 ss [i][j] = 1.0;
779 }
780 else
781 {
782 ss [i][j] = 0.0;
783 }
784 }
785 }
786
787 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
788 {
789 G4ThreeVector r = participants[i]->GetPosition();
790 rr[0][0] += r.y() * r.y() + r.z() * r.z();
791 rr[1][0] += - r.y() * r.x();
792 rr[2][0] += - r.z() * r.x();
793 rr[0][1] += - r.x() * r.y();
794 rr[1][1] += r.z() * r.z() + r.x() * r.x();
795 rr[2][1] += - r.z() * r.y();
796 rr[2][0] += - r.x() * r.z();
797 rr[2][1] += - r.y() * r.z();
798 rr[2][2] += r.x() * r.x() + r.y() * r.y();
799 }
800
801 for ( G4int i = 0 ; i < 3 ; i++ )
802 {
803 G4double x = rr [i][i];
804 for ( G4int j = 0 ; j < 3 ; j++ )
805 {
806 rr[i][j] = rr[i][j] / x;
807 ss[i][j] = ss[i][j] / x;
808 }
809 for ( G4int j = 0 ; j < 3 ; j++ )
810 {
811 if ( j != i )
812 {
813 G4double y = rr [j][i];
814 for ( G4int k = 0 ; k < 3 ; k++ )
815 {
816 rr[j][k] += -y * rr[i][k];
817 ss[j][k] += -y * ss[i][k];
818 }
819 }
820 }
821 }
822
823 G4double opl[3];
824 G4double rll[3];
825
826 rll[0] = ll.x();
827 rll[1] = ll.y();
828 rll[2] = ll.z();
829
830 for ( G4int i = 0 ; i < 3 ; i++ )
831 {
832 opl[i] = 0.0;
833
834 for ( G4int j = 0 ; j < 3 ; j++ )
835 {
836 opl[i] += ss[i][j]*rll[j];
837 }
838 }
839
840 for ( G4int i = 0 ; i < GetMassNumber() ; i++ )
841 {
842 G4ThreeVector p_i = participants[i]->GetMomentum() ;
843 G4ThreeVector ri = participants[i]->GetPosition() ;
844 G4ThreeVector opl_v ( opl[0] , opl[1] , opl[2] );
845
846 p_i += ri.cross(opl_v);
847
848 participants[i]->SetMomentum( p_i );
849 }
850
851}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
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 G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double diff2(const Hep3Vector &v) const
double y() const
Hep3Vector cross(const Hep3Vector &) const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4double GetBindingEnergy(const G4int A, const G4int Z)
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double A13(G4double A) const
Definition: G4Pow.cc:120
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4double GetTotalPotential()
G4double GetRHE(G4int i, G4int j)
void Cal2BodyQuantities()
G4ThreeVector GetFFr(G4int i)
G4double GetRHA(G4int i, G4int j)
void SetSystem(G4QMDSystem *aSystem)
G4double GetRR2(G4int i, G4int j)
G4ThreeVector GetFFp(G4int i)
G4int GetAtomicNumber()
Definition: G4QMDNucleus.cc:89
G4int GetMassNumber()
Definition: G4QMDNucleus.cc:68
G4double Get_clp()
G4double Get_cpc()
G4double Get_cdp()
G4double Get_epsx()
static G4QMDParameters * GetInstance()
G4double Get_hbc()
G4double Get_c0p()
G4double Get_gamm()
G4double Get_cpw()
G4double Get_c3p()
G4double Get_cph()
G4double Get_csp()
std::vector< G4QMDParticipant * > participants
Definition: G4QMDSystem.hh:72
void SetParticipant(G4QMDParticipant *particle)
Definition: G4QMDSystem.hh:51
#define DBL_MAX
Definition: templates.hh:62