hepmc - Diff between revs 96 and 103

Subversion Repositories:
Rev:
Only display areas with differences - Ignore whitespace
Rev 96 Rev 103
1 ////////////////////////////////////////////////////////////////////////// 1 //////////////////////////////////////////////////////////////////////////
2 // Matt.Dobbs@Cern.CH, September 1999 2 // Matt.Dobbs@Cern.CH, September 1999
3 // Updated: 7.1.2000 iterators complete and working! 3 // Updated: 7.1.2000 iterators complete and working!
4 // Updated: 10.1.2000 GenEvent::vertex, particle iterators are made 4 // Updated: 10.1.2000 GenEvent::vertex, particle iterators are made
5 //                    constant WRT this event ... note that 5 //                    constant WRT this event ... note that
6 //                    GenVertex::***_iterator is not const, since it must 6 //                    GenVertex::***_iterator is not const, since it must
7 //                    be able to return a mutable pointer to itself. 7 //                    be able to return a mutable pointer to itself.
8 // Updated: 08.2.2000 the event now holds a set of all attached vertices 8 // Updated: 08.2.2000 the event now holds a set of all attached vertices
9 //                    rather than just the roots of the graph 9 //                    rather than just the roots of the graph
10 // Event record for MC generators (for use at any stage of generation) 10 // Event record for MC generators (for use at any stage of generation)
11 ////////////////////////////////////////////////////////////////////////// 11 //////////////////////////////////////////////////////////////////////////
12 12
13 #include "HepMC/GenEvent.h" 13 #include "HepMC/GenEvent.h"
14 14
15 namespace HepMC { 15 namespace HepMC {
16 16
17     GenEvent::GenEvent( int signal_process_id, 17     GenEvent::GenEvent( int signal_process_id,
18                         int event_number, 18                         int event_number,
19                         GenVertex* signal_vertex, 19                         GenVertex* signal_vertex,
20                         const WeightContainer& weights, 20                         const WeightContainer& weights,
21                         const std::vector<long int>& random_states ) : 21                         const std::vector<long int>& random_states ) :
22         m_signal_process_id(signal_process_id), m_event_number(event_number), 22         m_signal_process_id(signal_process_id), m_event_number(event_number),
-   23         m_mpi(-1),
23         m_event_scale(-1), m_alphaQCD(-1), m_alphaQED(-1), 24         m_event_scale(-1), m_alphaQCD(-1), m_alphaQED(-1),
24         m_signal_process_vertex(signal_vertex), m_weights(weights), 25         m_signal_process_vertex(signal_vertex), m_weights(weights),
25         m_random_states(random_states), 26         m_random_states(random_states),
26         m_heavy_ion(0), 27         m_heavy_ion(0),
27         m_pdf_info(0) 28         m_pdf_info(0)
28     { 29     {
29         /// This constructor only allows null pointers to HeavyIon and PdfInfo 30         /// This constructor only allows null pointers to HeavyIon and PdfInfo
30         /// 31         ///
31         /// note: default values for m_event_scale, m_alphaQCD, m_alphaQED 32         /// note: default values for m_event_scale, m_alphaQCD, m_alphaQED
32         ///       are as suggested in hep-ph/0109068, "Generic Interface..." 33         ///       are as suggested in hep-ph/0109068, "Generic Interface..."
33         s_counter++; 34         s_counter++;
34     } 35     }
35 36
36     GenEvent::GenEvent( int signal_process_id, int event_number, 37     GenEvent::GenEvent( int signal_process_id, int event_number,
37                         GenVertex* signal_vertex, 38                         GenVertex* signal_vertex,
38                         const WeightContainer& weights, 39                         const WeightContainer& weights,
39                         const std::vector<long int>& random_states, 40                         const std::vector<long int>& random_states,
40                         const HeavyIon& ion, 41                         const HeavyIon& ion,
41                         const PdfInfo& pdf ) : 42                         const PdfInfo& pdf ) :
42         m_signal_process_id(signal_process_id), m_event_number(event_number), 43         m_signal_process_id(signal_process_id), m_event_number(event_number),
-   44         m_mpi(-1),
43         m_event_scale(-1), m_alphaQCD(-1), m_alphaQED(-1), 45         m_event_scale(-1), m_alphaQCD(-1), m_alphaQED(-1),
44         m_signal_process_vertex(signal_vertex), m_weights(weights), 46         m_signal_process_vertex(signal_vertex), m_weights(weights),
45         m_random_states(random_states), 47         m_random_states(random_states),
46         m_heavy_ion( new HeavyIon(ion) ), 48         m_heavy_ion( new HeavyIon(ion) ),
47         m_pdf_info( new PdfInfo(pdf) ) 49         m_pdf_info( new PdfInfo(pdf) )
48     { 50     {
49         /// GenEvent makes its own copy of HeavyIon and PdfInfo 51         /// GenEvent makes its own copy of HeavyIon and PdfInfo
50         /// 52         ///
51         /// note: default values for m_event_scale, m_alphaQCD, m_alphaQED 53         /// note: default values for m_event_scale, m_alphaQCD, m_alphaQED
52         ///       are as suggested in hep-ph/0109068, "Generic Interface..." 54         ///       are as suggested in hep-ph/0109068, "Generic Interface..."
53         s_counter++; 55         s_counter++;
54     } 56     }
55 57
56     GenEvent::GenEvent( const GenEvent& inevent ) 58     GenEvent::GenEvent( const GenEvent& inevent )
57     { 59     {
58         /// deep copy 60         /// deep copy
59         *this = inevent; 61         *this = inevent;
60         s_counter++; 62         s_counter++;
61     } 63     }
62 64
63     GenEvent::~GenEvent() 65     GenEvent::~GenEvent()
64     { 66     {
65         /// Deep destructor. 67         /// Deep destructor.
66         /// deletes all vertices/particles in this evt 68         /// deletes all vertices/particles in this evt
67         /// 69         ///
68         delete_all_vertices(); 70         delete_all_vertices();
69         delete m_heavy_ion; 71         delete m_heavy_ion;
70         delete m_pdf_info; 72         delete m_pdf_info;
71         s_counter--; 73         s_counter--;
72     } 74     }
73 75
74     GenEvent& GenEvent::operator=( const GenEvent& inevent ) 76     GenEvent& GenEvent::operator=( const GenEvent& inevent )
75     { 77     {
76         /// deep - makes a copy of all vertices! 78         /// deep - makes a copy of all vertices!
77         // 79         //
78         // 1. Delete all vertices attached to this 80         // 1. Delete all vertices attached to this
79         delete_all_vertices(); 81         delete_all_vertices();
80         //     82         //    
81         // 2. create a NEW copy of all vertices from inevent 83         // 2. create a NEW copy of all vertices from inevent
82         //    taking care to map new vertices onto the vertices being copied 84         //    taking care to map new vertices onto the vertices being copied
83         //    and add these new vertices to this event. 85         //    and add these new vertices to this event.
84         //    We do not use GenVertex::operator= because that would copy 86         //    We do not use GenVertex::operator= because that would copy
85         //    the attached particles as well. 87         //    the attached particles as well.
86         std::map<const GenVertex*,GenVertex*> map_in_to_new; 88         std::map<const GenVertex*,GenVertex*> map_in_to_new;
87         for ( GenEvent::vertex_const_iterator v = inevent.vertices_begin(); 89         for ( GenEvent::vertex_const_iterator v = inevent.vertices_begin();
88               v != inevent.vertices_end(); v++ ) { 90               v != inevent.vertices_end(); v++ ) {
89             GenVertex* newvertex = new GenVertex( 91             GenVertex* newvertex = new GenVertex(
90                 (*v)->position(), (*v)->id(), (*v)->weights() ); 92                 (*v)->position(), (*v)->id(), (*v)->weights() );
91             newvertex->suggest_barcode( (*v)->barcode() ); 93             newvertex->suggest_barcode( (*v)->barcode() );
92             map_in_to_new[*v] = newvertex; 94             map_in_to_new[*v] = newvertex;
93             add_vertex( newvertex ); 95             add_vertex( newvertex );
94         } 96         }
95         // 2.b copy the signal process vertex info. 97         // 2.b copy the signal process vertex info.
96         if ( inevent.signal_process_vertex() ) { 98         if ( inevent.signal_process_vertex() ) {
97             set_signal_process_vertex( 99             set_signal_process_vertex(
98                 map_in_to_new[inevent.signal_process_vertex()] ); 100                 map_in_to_new[inevent.signal_process_vertex()] );
99         } else set_signal_process_vertex( 0 ); 101         } else set_signal_process_vertex( 0 );
100         // 102         //
101         // 3. create a NEW copy of all particles from inevent 103         // 3. create a NEW copy of all particles from inevent
102         //    taking care to attach them to the appropriate 104         //    taking care to attach them to the appropriate
103         for ( GenEvent::particle_const_iterator p = inevent.particles_begin(); 105         for ( GenEvent::particle_const_iterator p = inevent.particles_begin();
104               p != inevent.particles_end(); p++ ) 106               p != inevent.particles_end(); p++ )
105         { 107         {
106             GenParticle* oldparticle = *p; 108             GenParticle* oldparticle = *p;
107             GenParticle* newparticle = new GenParticle(*oldparticle); 109             GenParticle* newparticle = new GenParticle(*oldparticle);
108             if ( oldparticle->end_vertex() ) { 110             if ( oldparticle->end_vertex() ) {
109                 map_in_to_new[ oldparticle->end_vertex() ]-> 111                 map_in_to_new[ oldparticle->end_vertex() ]->
110                                          add_particle_in(newparticle); 112                                          add_particle_in(newparticle);
111             } 113             }
112             if ( oldparticle->production_vertex() ) { 114             if ( oldparticle->production_vertex() ) {
113                 map_in_to_new[ oldparticle->production_vertex() ]-> 115                 map_in_to_new[ oldparticle->production_vertex() ]->
114                                          add_particle_out(newparticle); 116                                          add_particle_out(newparticle);
115             } 117             }
116         } 118         }
117         // 119         //
118         // 4. now that vtx/particles are copied, do everything else 120         // 4. now that vtx/particles are copied, do everything else
119         set_signal_process_id( inevent.signal_process_id() ); 121         set_signal_process_id( inevent.signal_process_id() );
120         set_event_number( inevent.event_number() ); 122         set_event_number( inevent.event_number() );
121         set_event_scale( inevent.event_scale() ); 123         set_event_scale( inevent.event_scale() );
122         set_alphaQCD( inevent.alphaQCD() ); 124         set_alphaQCD( inevent.alphaQCD() );
123         set_alphaQED( inevent.alphaQED() ); 125         set_alphaQED( inevent.alphaQED() );
124         set_random_states( inevent.random_states() ); 126         set_random_states( inevent.random_states() );
125         weights() = inevent.weights(); 127         weights() = inevent.weights();
126         // 128         //
127         // 5. copy these only if they are not null 129         // 5. copy these only if they are not null
128         if( inevent.heavy_ion() ) set_heavy_ion( *inevent.heavy_ion() ); 130         if( inevent.heavy_ion() ) set_heavy_ion( *inevent.heavy_ion() );
129         if( inevent.pdf_info() ) set_pdf_info( *inevent.pdf_info() ); 131         if( inevent.pdf_info() ) set_pdf_info( *inevent.pdf_info() );
130         return *this; 132         return *this;
131     } 133     }
132 134
133     void GenEvent::print( std::ostream& ostr ) const { 135     void GenEvent::print( std::ostream& ostr ) const {
134         /// dumps the content of this event to ostr 136         /// dumps the content of this event to ostr
135         ///   to dump to cout use: event.print(); 137         ///   to dump to cout use: event.print();
136         ///   if you want to write this event to file outfile.txt you could use: 138         ///   if you want to write this event to file outfile.txt you could use:
137         ///      std::ofstream outfile("outfile.txt"); event.print( outfile ); 139         ///      std::ofstream outfile("outfile.txt"); event.print( outfile );
138         ostr << "________________________________________" 140         ostr << "________________________________________"
139              << "________________________________________\n"; 141              << "________________________________________\n";
140         ostr << "GenEvent: #" << event_number() 142         ostr << "GenEvent: #" << event_number()
141              << " ID=" << signal_process_id() 143              << " ID=" << signal_process_id()
142              << " SignalProcessGenVertex Barcode: " 144              << " SignalProcessGenVertex Barcode: "
143              << ( signal_process_vertex() ? signal_process_vertex()->barcode() 145              << ( signal_process_vertex() ? signal_process_vertex()->barcode()
144                   : 0 ) 146                   : 0 )
145              << "\n"; 147              << "\n";
146         ostr << " Current Memory Usage: " 148         ostr << " Current Memory Usage: "
147              << GenEvent::counter() << " events, " 149              << GenEvent::counter() << " events, "
148              << GenVertex::counter() << " vertices, " 150              << GenVertex::counter() << " vertices, "
149              << GenParticle::counter() << " particles.\n"; 151              << GenParticle::counter() << " particles.\n";
150         ostr << " Entries this event: " << vertices_size() << " vertices, " 152         ostr << " Entries this event: " << vertices_size() << " vertices, "
151              << particles_size() << " particles.\n"; 153              << particles_size() << " particles.\n";
152         // Random State 154         // Random State
153         ostr << " RndmState(" << m_random_states.size() << ")="; 155         ostr << " RndmState(" << m_random_states.size() << ")=";
154         for ( std::vector<long int>::const_iterator rs 156         for ( std::vector<long int>::const_iterator rs
155                   = m_random_states.begin(); 157                   = m_random_states.begin();
156               rs != m_random_states.end(); rs++ ) { ostr << *rs << " "; } 158               rs != m_random_states.end(); rs++ ) { ostr << *rs << " "; }
157         ostr << "\n"; 159         ostr << "\n";
158         // Weights 160         // Weights
159         ostr << " Wgts(" << weights().size() << ")="; 161         ostr << " Wgts(" << weights().size() << ")=";
160         for ( WeightContainer::const_iterator wgt = weights().begin(); 162         for ( WeightContainer::const_iterator wgt = weights().begin();
161               wgt != weights().end(); wgt++ ) { ostr << *wgt << " "; } 163               wgt != weights().end(); wgt++ ) { ostr << *wgt << " "; }
162         ostr << "\n"; 164         ostr << "\n";
163         ostr << " EventScale " << event_scale() 165         ostr << " EventScale " << event_scale()
164              << " [energy] \t alphaQCD=" << alphaQCD() 166              << " [energy] \t alphaQCD=" << alphaQCD()
165              << "\t alphaQED=" << alphaQED() << std::endl; 167              << "\t alphaQED=" << alphaQED() << std::endl;
166         // print a legend to describe the particle info 168         // print a legend to describe the particle info
167         ostr << "                                    GenParticle Legend\n"; 169         ostr << "                                    GenParticle Legend\n";
168         ostr  << "        Barcode   PDG ID      " 170         ostr  << "        Barcode   PDG ID      "
169               << "( Px,       Py,       Pz,     E )" 171               << "( Px,       Py,       Pz,     E )"
170               << " Stat  DecayVtx\n"; 172               << " Stat  DecayVtx\n";
171         ostr << "________________________________________" 173         ostr << "________________________________________"
172              << "________________________________________\n"; 174              << "________________________________________\n";
173         // Print all Vertices 175         // Print all Vertices
174         for ( GenEvent::vertex_const_iterator vtx = this->vertices_begin(); 176         for ( GenEvent::vertex_const_iterator vtx = this->vertices_begin();
175               vtx != this->vertices_end(); ++vtx ) { 177               vtx != this->vertices_end(); ++vtx ) {
176             (*vtx)->print(ostr); 178             (*vtx)->print(ostr);
177         } 179         }
178         ostr << "________________________________________" 180         ostr << "________________________________________"
179              << "________________________________________" << std::endl; 181              << "________________________________________" << std::endl;
180     } 182     }
181 183
182     bool GenEvent::add_vertex( GenVertex* vtx ) { 184     bool GenEvent::add_vertex( GenVertex* vtx ) {
183         /// returns true if successful - generally will only return false 185         /// returns true if successful - generally will only return false
184         /// if the inserted vertex is already included in the event. 186         /// if the inserted vertex is already included in the event.
185         if ( !vtx ) return 0; 187         if ( !vtx ) return 0;
186         // if vtx previously pointed to another GenEvent, remove it from that 188         // if vtx previously pointed to another GenEvent, remove it from that
187         // GenEvent's list 189         // GenEvent's list
188         if ( vtx->parent_event() && vtx->parent_event() != this ) { 190         if ( vtx->parent_event() && vtx->parent_event() != this ) {
189             bool remove_status = vtx->parent_event()->remove_vertex( vtx ); 191             bool remove_status = vtx->parent_event()->remove_vertex( vtx );
190             if ( !remove_status ) {             192             if ( !remove_status ) {            
191                 std::cerr << "GenEvent::add_vertex ERROR " 193                 std::cerr << "GenEvent::add_vertex ERROR "
192                           << "GenVertex::parent_event points to \n" 194                           << "GenVertex::parent_event points to \n"
193                           << "an event that does not point back to the " 195                           << "an event that does not point back to the "
194                           << "GenVertex. \n This probably indicates a deeper " 196                           << "GenVertex. \n This probably indicates a deeper "
195                           << "problem. " << std::endl; 197                           << "problem. " << std::endl;
196             } 198             }
197         } 199         }
198         // 200         //
199         // setting the vertex parent also inserts the vertex into this 201         // setting the vertex parent also inserts the vertex into this
200         // event 202         // event
201         vtx->set_parent_event_( this ); 203         vtx->set_parent_event_( this );
202         return ( m_vertex_barcodes.count(vtx->barcode()) ? true : false ); 204         return ( m_vertex_barcodes.count(vtx->barcode()) ? true : false );
203     } 205     }
204 206
205     bool GenEvent::remove_vertex( GenVertex* vtx ) { 207     bool GenEvent::remove_vertex( GenVertex* vtx ) {
206         /// this removes vtx from the event but does NOT delete it. 208         /// this removes vtx from the event but does NOT delete it.
207         /// returns True if an entry vtx existed in the table and was erased 209         /// returns True if an entry vtx existed in the table and was erased
208         if ( m_signal_process_vertex == vtx ) m_signal_process_vertex = 0; 210         if ( m_signal_process_vertex == vtx ) m_signal_process_vertex = 0;
209         if ( vtx->parent_event() == this ) vtx->set_parent_event_( 0 ); 211         if ( vtx->parent_event() == this ) vtx->set_parent_event_( 0 );
210         return ( m_vertex_barcodes.count(vtx->barcode()) ? false : true ); 212         return ( m_vertex_barcodes.count(vtx->barcode()) ? false : true );
211     } 213     }
212 214
213     void GenEvent::clear() 215     void GenEvent::clear()
214     { 216     {
215         /// remove all information from the event 217         /// remove all information from the event
216         /// deletes all vertices/particles in this evt 218         /// deletes all vertices/particles in this evt
217         /// 219         ///
218         delete_all_vertices(); 220         delete_all_vertices();
219         delete m_heavy_ion; 221         delete m_heavy_ion;
220         delete m_pdf_info; 222         delete m_pdf_info;
221         m_signal_process_id = 0; 223         m_signal_process_id = 0;
222         m_event_number = 0; 224         m_event_number = 0;
223         m_event_scale = -1; 225         m_event_scale = -1;
224         m_alphaQCD = -1; 226         m_alphaQCD = -1;
225         m_alphaQED = -1; 227         m_alphaQED = -1;
226         m_weights = std::vector<double>(); 228         m_weights = std::vector<double>();
227         m_random_states = std::vector<long int>(); 229         m_random_states = std::vector<long int>();
228         // error check just to be safe 230         // error check just to be safe
229         if ( m_vertex_barcodes.size() != 0 231         if ( m_vertex_barcodes.size() != 0
230              || m_particle_barcodes.size() != 0 ) { 232              || m_particle_barcodes.size() != 0 ) {
231             std::cerr << "GenEvent::clear() strange result ... \n" 233             std::cerr << "GenEvent::clear() strange result ... \n"
232                       << "either the particle and/or the vertex map isn't empty" << std::endl; 234                       << "either the particle and/or the vertex map isn't empty" << std::endl;
233             std::cerr << "Number vtx,particle the event after deleting = " 235             std::cerr << "Number vtx,particle the event after deleting = "
234                       << m_vertex_barcodes.size() << "  " 236                       << m_vertex_barcodes.size() << "  "
235                       << m_particle_barcodes.size() << std::endl; 237                       << m_particle_barcodes.size() << std::endl;
236             std::cerr << "Total Number vtx,particle in memory " 238             std::cerr << "Total Number vtx,particle in memory "
237                       << "after method called = " 239                       << "after method called = "
238                       << GenVertex::counter() << "\t" 240                       << GenVertex::counter() << "\t"
239                       << GenParticle::counter() << std::endl; 241                       << GenParticle::counter() << std::endl;
240         } 242         }
241         return; 243         return;
242     } 244     }
243     245    
244     void GenEvent::delete_all_vertices() { 246     void GenEvent::delete_all_vertices() {
245         /// deletes all vertices in the vertex container 247         /// deletes all vertices in the vertex container
246         /// (i.e. all vertices owned by this event) 248         /// (i.e. all vertices owned by this event)
247         /// The vertices are the "owners" of the particles, so as we delete 249         /// The vertices are the "owners" of the particles, so as we delete
248         ///   the vertices, the vertex desctructors are automatically 250         ///   the vertices, the vertex desctructors are automatically
249         ///   deleting their particles. 251         ///   deleting their particles.
250         if ( vertices_empty() ) return; 252         if ( vertices_empty() ) return;
251         // delete each vertex individually (this deletes particles as well) 253         // delete each vertex individually (this deletes particles as well)
252         while ( m_vertex_barcodes.begin() != m_vertex_barcodes.end() ) { 254         while ( m_vertex_barcodes.begin() != m_vertex_barcodes.end() ) {
253             GenVertex* vtx = ( m_vertex_barcodes.begin() )->second; 255             GenVertex* vtx = ( m_vertex_barcodes.begin() )->second;
254             m_vertex_barcodes.erase( m_vertex_barcodes.begin() ); 256             m_vertex_barcodes.erase( m_vertex_barcodes.begin() );
255             delete vtx; 257             delete vtx;
256         } 258         }
257         // 259         //
258         // Error checking: 260         // Error checking:
259         if ( !vertices_empty() || ! particles_empty() ) { 261         if ( !vertices_empty() || ! particles_empty() ) {
260             std::cerr << "GenEvent::delete_all_vertices strange result ... " 262             std::cerr << "GenEvent::delete_all_vertices strange result ... "
261                       << "after deleting all vertices, \nthe particle and " 263                       << "after deleting all vertices, \nthe particle and "
262                       << "vertex maps aren't empty.\n  This probably " 264                       << "vertex maps aren't empty.\n  This probably "
263                       << "indicates deeper problems or memory leak in the " 265                       << "indicates deeper problems or memory leak in the "
264                       << "code." << std::endl; 266                       << "code." << std::endl;
265             std::cerr << "Number vtx,particle the event after deleting = " 267             std::cerr << "Number vtx,particle the event after deleting = "
266                       << m_vertex_barcodes.size() << "  " 268                       << m_vertex_barcodes.size() << "  "
267                       << m_particle_barcodes.size() << std::endl; 269                       << m_particle_barcodes.size() << std::endl;
268             std::cerr << "Total Number vtx,particle in memory " 270             std::cerr << "Total Number vtx,particle in memory "
269                       << "after method called = " 271                       << "after method called = "
270                       << GenVertex::counter() << "\t" 272                       << GenVertex::counter() << "\t"
271                       << GenParticle::counter() << std::endl; 273                       << GenParticle::counter() << std::endl;
272         } 274         }
273     } 275     }
274     276    
275     bool GenEvent::set_barcode( GenParticle* p, int suggested_barcode ) 277     bool GenEvent::set_barcode( GenParticle* p, int suggested_barcode )
276     { 278     {
277         if ( p->parent_event() != this ) { 279         if ( p->parent_event() != this ) {
278             std::cerr << "GenEvent::set_barcode attempted, but the argument's" 280             std::cerr << "GenEvent::set_barcode attempted, but the argument's"
279                       << "\n parent_event is not this ... request rejected." 281                       << "\n parent_event is not this ... request rejected."
280                       << std::endl; 282                       << std::endl;
281             return false; 283             return false;
282         } 284         }
283         // M.Dobbs  Nov 4, 2002 285         // M.Dobbs  Nov 4, 2002
284         // First we must check to see if the particle already has a 286         // First we must check to see if the particle already has a
285         // barcode which is different from the suggestion. If yes, we 287         // barcode which is different from the suggestion. If yes, we
286         // remove it from the particle map. 288         // remove it from the particle map.
287         if ( p->barcode() != 0 && p->barcode() != suggested_barcode ) { 289         if ( p->barcode() != 0 && p->barcode() != suggested_barcode ) {
288             if ( m_particle_barcodes.count(p->barcode()) && 290             if ( m_particle_barcodes.count(p->barcode()) &&
289                  m_particle_barcodes[p->barcode()] == p ) { 291                  m_particle_barcodes[p->barcode()] == p ) {
290                 m_particle_barcodes.erase( p->barcode() ); 292                 m_particle_barcodes.erase( p->barcode() );
291             } 293             }
292             // At this point either the particle is NOT in 294             // At this point either the particle is NOT in
293             // m_particle_barcodes, or else it is in the map, but 295             // m_particle_barcodes, or else it is in the map, but
294             // already with the suggested barcode. 296             // already with the suggested barcode.
295         } 297         }
296         // 298         //
297         // First case --- a valid barcode has been suggested 299         // First case --- a valid barcode has been suggested
298         //     (valid barcodes are numbers greater than zero) 300         //     (valid barcodes are numbers greater than zero)
299         bool insert_success = true; 301         bool insert_success = true;
300         if ( suggested_barcode > 0 ) { 302         if ( suggested_barcode > 0 ) {
301             if ( m_particle_barcodes.count(suggested_barcode) ) { 303             if ( m_particle_barcodes.count(suggested_barcode) ) {
302                 // the suggested_barcode is already used. 304                 // the suggested_barcode is already used.
303                 if ( m_particle_barcodes[suggested_barcode] == p ) { 305                 if ( m_particle_barcodes[suggested_barcode] == p ) {
304                     // but it was used for this particle ... so everythings ok 306                     // but it was used for this particle ... so everythings ok
305                     p->set_barcode_( suggested_barcode ); 307                     p->set_barcode_( suggested_barcode );
306                     return true; 308                     return true;
307                 } 309                 }
308                 insert_success = false; 310                 insert_success = false;
309                 suggested_barcode = 0; 311                 suggested_barcode = 0;
310             } else { // suggested barcode is OK, proceed to insert 312             } else { // suggested barcode is OK, proceed to insert
311                 m_particle_barcodes[suggested_barcode] = p; 313                 m_particle_barcodes[suggested_barcode] = p;
312                 p->set_barcode_( suggested_barcode ); 314                 p->set_barcode_( suggested_barcode );
313                 return true; 315                 return true;
314             } 316             }
315         } 317         }
316         // 318         //
317         // Other possibility -- a valid barcode has not been suggested, 319         // Other possibility -- a valid barcode has not been suggested,
318         //    so one is automatically generated 320         //    so one is automatically generated
319         if ( suggested_barcode < 0 ) insert_success = false; 321         if ( suggested_barcode < 0 ) insert_success = false;
320         if ( suggested_barcode <= 0 ) { 322         if ( suggested_barcode <= 0 ) {
321             if ( !m_particle_barcodes.empty() ) { 323             if ( !m_particle_barcodes.empty() ) {
322                 // in this case we find the highest barcode that was used, 324                 // in this case we find the highest barcode that was used,
323                 // and increment it by 1 325                 // and increment it by 1
324                 suggested_barcode = m_particle_barcodes.rbegin()->first; 326                 suggested_barcode = m_particle_barcodes.rbegin()->first;
325                 ++suggested_barcode; 327                 ++suggested_barcode;
326             } 328             }
327             // For the automatically assigned barcodes, the first one 329             // For the automatically assigned barcodes, the first one
328             //   we use is 10001 ... this is just because when we read 330             //   we use is 10001 ... this is just because when we read
329             //   events from HEPEVT, we will suggest their index as the 331             //   events from HEPEVT, we will suggest their index as the
330             //   barcode, and that index has maximum value 10000. 332             //   barcode, and that index has maximum value 10000.
331             //  This way we will immediately be able to recognize the hepevt 333             //  This way we will immediately be able to recognize the hepevt
332             //   particles from the auto-assigned ones. 334             //   particles from the auto-assigned ones.
333             if ( suggested_barcode <= 10000 ) suggested_barcode = 10001; 335             if ( suggested_barcode <= 10000 ) suggested_barcode = 10001;
334         } 336         }
335         // At this point we should have a valid barcode 337         // At this point we should have a valid barcode
336         if ( m_particle_barcodes.count(suggested_barcode) ) { 338         if ( m_particle_barcodes.count(suggested_barcode) ) {
337             std::cerr << "GenEvent::set_barcode ERROR, this should never " 339             std::cerr << "GenEvent::set_barcode ERROR, this should never "
338                       << "happen \n report bug to matt.dobbs@cern.ch" 340                       << "happen \n report bug to matt.dobbs@cern.ch"
339                       << std::endl; 341                       << std::endl;
340         } 342         }
341         m_particle_barcodes[suggested_barcode] = p; 343         m_particle_barcodes[suggested_barcode] = p;
342         p->set_barcode_( suggested_barcode ); 344         p->set_barcode_( suggested_barcode );
343         return insert_success; 345         return insert_success;
344     } 346     }
345 347
346     bool GenEvent::set_barcode( GenVertex* v, int suggested_barcode ) 348     bool GenEvent::set_barcode( GenVertex* v, int suggested_barcode )
347     { 349     {
348         if ( v->parent_event() != this ) { 350         if ( v->parent_event() != this ) {
349             std::cerr << "GenEvent::set_barcode attempted, but the argument's" 351             std::cerr << "GenEvent::set_barcode attempted, but the argument's"
350                       << "\n parent_event is not this ... request rejected." 352                       << "\n parent_event is not this ... request rejected."
351                       << std::endl; 353                       << std::endl;
352             return false; 354             return false;
353         } 355         }
354         // M.Dobbs Nov 4, 2002 356         // M.Dobbs Nov 4, 2002
355         // First we must check to see if the vertex already has a 357         // First we must check to see if the vertex already has a
356         // barcode which is different from the suggestion. If yes, we 358         // barcode which is different from the suggestion. If yes, we
357         // remove it from the vertex map. 359         // remove it from the vertex map.
358         if ( v->barcode() != 0 && v->barcode() != suggested_barcode ) { 360         if ( v->barcode() != 0 && v->barcode() != suggested_barcode ) {
359             if ( m_vertex_barcodes.count(v->barcode()) && 361             if ( m_vertex_barcodes.count(v->barcode()) &&
360                  m_vertex_barcodes[v->barcode()] == v ) { 362                  m_vertex_barcodes[v->barcode()] == v ) {
361                 m_vertex_barcodes.erase( v->barcode() ); 363                 m_vertex_barcodes.erase( v->barcode() );
362             } 364             }
363             // At this point either the vertex is NOT in 365             // At this point either the vertex is NOT in
364             // m_vertex_barcodes, or else it is in the map, but 366             // m_vertex_barcodes, or else it is in the map, but
365             // already with the suggested barcode. 367             // already with the suggested barcode.
366         } 368         }
367         369        
368         // 370         //
369         // First case --- a valid barcode has been suggested 371         // First case --- a valid barcode has been suggested
370         //     (valid barcodes are numbers greater than zero) 372         //     (valid barcodes are numbers greater than zero)
371         bool insert_success = true; 373         bool insert_success = true;
372         if ( suggested_barcode < 0 ) { 374         if ( suggested_barcode < 0 ) {
373             if ( m_vertex_barcodes.count(suggested_barcode) ) { 375             if ( m_vertex_barcodes.count(suggested_barcode) ) {
374                 // the suggested_barcode is already used. 376                 // the suggested_barcode is already used.
375                 if ( m_vertex_barcodes[suggested_barcode] == v ) { 377                 if ( m_vertex_barcodes[suggested_barcode] == v ) {
376                     // but it was used for this vertex ... so everythings ok 378                     // but it was used for this vertex ... so everythings ok
377                     v->set_barcode_( suggested_barcode ); 379                     v->set_barcode_( suggested_barcode );
378                     return true; 380                     return true;
379                 } 381                 }
380                 insert_success = false; 382                 insert_success = false;
381                 suggested_barcode = 0; 383                 suggested_barcode = 0;
382             } else { // suggested barcode is OK, proceed to insert 384             } else { // suggested barcode is OK, proceed to insert
383                 m_vertex_barcodes[suggested_barcode] = v; 385                 m_vertex_barcodes[suggested_barcode] = v;
384                 v->set_barcode_( suggested_barcode ); 386                 v->set_barcode_( suggested_barcode );
385                 return true; 387                 return true;
386             } 388             }
387         } 389         }
388         // 390         //
389         // Other possibility -- a valid barcode has not been suggested, 391         // Other possibility -- a valid barcode has not been suggested,
390         //    so one is automatically generated 392         //    so one is automatically generated
391         if ( suggested_barcode > 0 ) insert_success = false; 393         if ( suggested_barcode > 0 ) insert_success = false;
392         if ( suggested_barcode >= 0 ) { 394         if ( suggested_barcode >= 0 ) {
393             if ( !m_vertex_barcodes.empty() ) { 395             if ( !m_vertex_barcodes.empty() ) {
394                 // in this case we find the highest barcode that was used, 396                 // in this case we find the highest barcode that was used,
395                 // and increment it by 1, (vertex barcodes are negative) 397                 // and increment it by 1, (vertex barcodes are negative)
396                 suggested_barcode = m_vertex_barcodes.rbegin()->first; 398                 suggested_barcode = m_vertex_barcodes.rbegin()->first;
397                 --suggested_barcode; 399                 --suggested_barcode;
398             } 400             }
399             if ( suggested_barcode >= 0 ) suggested_barcode = -1; 401             if ( suggested_barcode >= 0 ) suggested_barcode = -1;
400         } 402         }
401         // At this point we should have a valid barcode 403         // At this point we should have a valid barcode
402         if ( m_vertex_barcodes.count(suggested_barcode) ) { 404         if ( m_vertex_barcodes.count(suggested_barcode) ) {
403             std::cerr << "GenEvent::set_barcode ERROR, this should never " 405             std::cerr << "GenEvent::set_barcode ERROR, this should never "
404                       << "happen \n report bug to matt.dobbs@cern.ch" 406                       << "happen \n report bug to matt.dobbs@cern.ch"
405                       << std::endl; 407                       << std::endl;
406         } 408         }
407         m_vertex_barcodes[suggested_barcode] = v; 409         m_vertex_barcodes[suggested_barcode] = v;
408         v->set_barcode_( suggested_barcode ); 410         v->set_barcode_( suggested_barcode );
409         return insert_success; 411         return insert_success;
410     } 412     }
411 413
412     ///////////// 414     /////////////
413     // Static  // 415     // Static  //
414     ///////////// 416     /////////////
415     unsigned int GenEvent::counter() { return s_counter; } 417     unsigned int GenEvent::counter() { return s_counter; }
416     unsigned int GenEvent::s_counter = 0; 418     unsigned int GenEvent::s_counter = 0;
417 419
418 } // HepMC 420 } // HepMC
419 421
420 422
421 423
422   424