source: trunk/source/processes/hadronic/models/qmd/src/G4QMDSystem.cc@ 1201

Last change on this file since 1201 was 962, checked in by garnier, 17 years ago

update processes

File size: 3.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// 081120 Add InsertParticipant Method by T. Koi
27
28#include "G4QMDSystem.hh"
29#include <iomanip>
30
31G4QMDSystem::G4QMDSystem()
32{
33 participants.clear();
34 numberOfCollision = 0;
35}
36
37
38
39G4QMDSystem::~G4QMDSystem()
40{
41 this->Clear();
42}
43
44
45// Insert nucleus to current system;
46void G4QMDSystem::SetSystem ( G4QMDSystem* nucleus , G4ThreeVector dp , G4ThreeVector dr )
47{
48 std::vector< G4QMDParticipant* >::iterator it;
49 for ( it = nucleus->participants.begin() ; it != nucleus->participants.end() ; it++ )
50 {
51 G4ThreeVector r = (*it)->GetPosition() + dr;
52 (*it)->SetPosition ( r );
53 G4ThreeVector p = (*it)->GetMomentum() + dp;
54 (*it)->SetMomentum ( p );
55 this->SetParticipant( *it );
56 }
57}
58
59void G4QMDSystem::SubtractSystem ( G4QMDSystem* nucleus )
60{
61
62 for ( G4int i = 0 ; i < nucleus->GetTotalNumberOfParticipant() ; i++ )
63 {
64 participants.erase ( std::find ( participants.begin() , participants.end() , nucleus->GetParticipant( i ) ) );
65 }
66}
67
68void G4QMDSystem::Clear ()
69{
70 for ( G4int i = 0 ; i < this->GetTotalNumberOfParticipant() ; i++ )
71 {
72 delete participants[i];
73 }
74 participants.clear();
75}
76
77
78
79void G4QMDSystem::ShowParticipants()
80{
81 G4ThreeVector p_sum( 0.0 );
82 std::vector< G4QMDParticipant* >::iterator it;
83 G4cout << "Momentum and Position of each participant " << G4endl;
84 G4int i = 0;
85 for ( it = participants.begin() ; it != participants.end() ; it++ )
86 {
87 G4cout << i
88 << " "
89 << (*it)->GetDefinition()->GetParticleName()
90 << " "
91 << std::setprecision( 8 )
92 << (*it)->GetMomentum()
93 << " "
94 << (*it)->GetPosition()
95 << G4endl;
96 p_sum += (*it)->GetMomentum();
97 i++;
98 }
99 G4cout << "Sum upped Momentum and mag " << p_sum << " " << p_sum.mag() << G4endl;
100}
101
102
103
104void G4QMDSystem::InsertParticipant ( G4QMDParticipant* particle , G4int n )
105{
106
107 if ( (size_t) n > participants.size()+1 )
108 G4cout << "G4QMDSystem::InsertParticipant size error" << G4endl;
109
110 std::vector< G4QMDParticipant* >::iterator it;
111 it = participants.begin();
112
113 for ( G4int i = 0; i < n ; i++ )
114 it++;
115
116 participants.insert( it, particle );
117}
Note: See TracBrowser for help on using the repository browser.