Fork me on GitHub

source: svn/trunk/Utilities/ExRootAnalysis/interface/LHEF.h@ 769

Last change on this file since 769 was 3, checked in by Xavier Rouby, 16 years ago

first commit

File size: 23.0 KB
Line 
1// -*- C++ -*-
2#ifndef THEPEG_LHEF_H
3#define THEPEG_LHEF_H
4//
5// This is the declaration of the Les Houches Event File classes.
6//
7
8
9#include <iostream>
10#include <iomanip>
11#include <sstream>
12#include <fstream>
13#include <string>
14#include <vector>
15#include <utility>
16#include <stdexcept>
17
18namespace LHEF {
19
20/**
21 * The HEPRUP class is a simple container corresponding to the Les Houches
22 * accord (<A HREF="http://arxiv.org/abs/hep-ph/0109068">hep-ph/0109068</A>)
23 * common block with the same name. The members are named in the same
24 * way as in the common block. However, fortran arrays are represented
25 * by vectors, except for the arrays of length two which are
26 * represented by pair objects.
27 */
28class HEPRUP {
29
30public:
31
32 /** @name Standard constructors and destructors. */
33 //@{
34 /**
35 * Default constructor.
36 */
37 HEPRUP()
38 : IDWTUP(0), NPRUP(0) {}
39
40 /**
41 * Copy-constructor.
42 */
43 HEPRUP(const HEPRUP & x)
44 : IDBMUP(x.IDBMUP), EBMUP(x.EBMUP),
45 PDFGUP(x.PDFGUP), PDFSUP(x.PDFSUP), IDWTUP(x.IDWTUP),
46 NPRUP(x.NPRUP), XSECUP(x.XSECUP), XERRUP(x.XERRUP),
47 XMAXUP(x.XMAXUP), LPRUP(x.LPRUP) {}
48
49
50 /**
51 * Assignment operator.
52 */
53 HEPRUP & operator=(const HEPRUP & x) {
54 IDBMUP = x.IDBMUP;
55 EBMUP = x.EBMUP;
56 PDFGUP = x.PDFGUP;
57 PDFSUP = x.PDFSUP;
58 IDWTUP = x.IDWTUP;
59 NPRUP = x.NPRUP;
60 XSECUP = x.XSECUP;
61 XERRUP = x.XERRUP;
62 XMAXUP = x.XMAXUP;
63 LPRUP = x.LPRUP;
64 return *this;
65 }
66
67 /**
68 * Destructor.
69 */
70 ~HEPRUP() {}
71 //@}
72
73public:
74
75 /**
76 * Set the NPRUP variable, corresponding to the number of
77 * sub-processes, to \a nrup, and resize all relevant vectors
78 * accordingly.
79 */
80 void resize(int nrup) {
81 NPRUP = nrup;
82 resize();
83 }
84
85 /**
86 * Assuming the NPRUP variable, corresponding to the number of
87 * sub-processes, is correctly set, resize the relevant vectors
88 * accordingly.
89 */
90 void resize() {
91 XSECUP.resize(NPRUP);
92 XERRUP.resize(NPRUP);
93 XMAXUP.resize(NPRUP);
94 LPRUP.resize(NPRUP);
95 }
96
97public:
98
99 /**
100 * PDG id's of beam particles. (first/second is in +/-z direction).
101 */
102 std::pair<long,long> IDBMUP;
103
104 /**
105 * Energy of beam particles given in GeV.
106 */
107 std::pair<double,double> EBMUP;
108
109 /**
110 * The author group for the PDF used for the beams according to the
111 * PDFLib specification.
112 */
113 std::pair<int,int> PDFGUP;
114
115 /**
116 * The id number the PDF used for the beams according to the
117 * PDFLib specification.
118 */
119 std::pair<int,int> PDFSUP;
120
121 /**
122 * Master switch indicating how the ME generator envisages the
123 * events weights should be interpreted according to the Les Houches
124 * accord.
125 */
126 int IDWTUP;
127
128 /**
129 * The number of different subprocesses in this file.
130 */
131 int NPRUP;
132
133 /**
134 * The cross sections for the different subprocesses in pb.
135 */
136 std::vector<double> XSECUP;
137
138 /**
139 * The statistical error in the cross sections for the different
140 * subprocesses in pb.
141 */
142 std::vector<double> XERRUP;
143
144 /**
145 * The maximum event weights (in HEPEUP::XWGTUP) for different
146 * subprocesses.
147 */
148 std::vector<double> XMAXUP;
149
150 /**
151 * The subprocess code for the different subprocesses.
152 */
153 std::vector<int> LPRUP;
154
155};
156
157
158/**
159 * The HEPEUP class is a simple container corresponding to the Les Houches accord
160 * (<A HREF="http://arxiv.org/abs/hep-ph/0109068">hep-ph/0109068</A>)
161 * common block with the same name. The members are named in the same
162 * way as in the common block. However, fortran arrays are represented
163 * by vectors, except for the arrays of length two which are
164 * represented by pair objects.
165 */
166class HEPEUP {
167
168public:
169
170 /** @name Standard constructors and destructors. */
171 //@{
172 /**
173 * Default constructor.
174 */
175 HEPEUP()
176 : NUP(0), IDPRUP(0), XWGTUP(0.0), XPDWUP(0.0, 0.0),
177 SCALUP(0.0), AQEDUP(0.0), AQCDUP(0.0) {}
178
179 /**
180 * Copy-constructor.
181 */
182 HEPEUP(const HEPEUP & x)
183 : NUP(x.NUP), IDPRUP(x.IDPRUP), XWGTUP(x.XWGTUP), XPDWUP(x.XPDWUP),
184 SCALUP(x.SCALUP), AQEDUP(x.AQEDUP), AQCDUP(x.AQCDUP), IDUP(x.IDUP),
185 ISTUP(x.ISTUP), MOTHUP(x.MOTHUP), ICOLUP(x.ICOLUP),
186 PUP(x.PUP), VTIMUP(x.VTIMUP), SPINUP(x.SPINUP) {}
187
188 /**
189 * Assignment operator.
190 */
191 HEPEUP & operator=(const HEPEUP & x) {
192 NUP = x.NUP;
193 IDPRUP = x.IDPRUP;
194 XWGTUP = x.XWGTUP;
195 XPDWUP = x.XPDWUP;
196 SCALUP = x.SCALUP;
197 AQEDUP = x.AQEDUP;
198 AQCDUP = x.AQCDUP;
199 IDUP = x.IDUP;
200 ISTUP = x.ISTUP;
201 MOTHUP = x.MOTHUP;
202 ICOLUP = x.ICOLUP;
203 PUP = x.PUP;
204 VTIMUP = x.VTIMUP;
205 SPINUP = x.SPINUP;
206 return *this;
207 }
208
209
210 /**
211 * Destructor.
212 */
213 ~HEPEUP() {};
214 //@}
215
216public:
217
218 /**
219 * Set the NUP variable, corresponding to the number of particles in
220 * the current event, to \a nup, and resize all relevant vectors
221 * accordingly.
222 */
223 void resize(int nup) {
224 NUP = nup;
225 resize();
226 }
227
228 /**
229 * Assuming the NUP variable, corresponding to the number of
230 * particles in the current event, is correctly set, resize the
231 * relevant vectors accordingly.
232 */
233 void resize() {
234 IDUP.resize(NUP);
235 ISTUP.resize(NUP);
236 MOTHUP.resize(NUP);
237 ICOLUP.resize(NUP);
238 PUP.resize(NUP, std::vector<double>(5));
239 VTIMUP.resize(NUP);
240 SPINUP.resize(NUP);
241 }
242
243public:
244
245 /**
246 * The number of particle entries in the current event.
247 */
248 int NUP;
249
250 /**
251 * The subprocess code for this event (as given in LPRUP).
252 */
253 int IDPRUP;
254
255 /**
256 * The weight for this event.
257 */
258 double XWGTUP;
259
260 /**
261 * The PDF weights for the two incoming partons. Note that this
262 * variable is not present in the current LesHouches accord
263 * (<A HREF="http://arxiv.org/abs/hep-ph/0109068">hep-ph/0109068</A>),
264 * hopefully it will be present in a future accord.
265 */
266 std::pair<double,double> XPDWUP;
267
268 /**
269 * The scale in GeV used in the calculation of the PDF's in this
270 * event.
271 */
272 double SCALUP;
273
274 /**
275 * The value of the QED coupling used in this event.
276 */
277 double AQEDUP;
278
279 /**
280 * The value of the QCD coupling used in this event.
281 */
282 double AQCDUP;
283
284 /**
285 * The PDG id's for the particle entries in this event.
286 */
287 std::vector<long> IDUP;
288
289 /**
290 * The status codes for the particle entries in this event.
291 */
292 std::vector<int> ISTUP;
293
294 /**
295 * Indices for the first and last mother for the particle entries in
296 * this event.
297 */
298 std::vector< std::pair<int,int> > MOTHUP;
299
300 /**
301 * The colour-line indices (first(second) is (anti)colour) for the
302 * particle entries in this event.
303 */
304 std::vector< std::pair<int,int> > ICOLUP;
305
306 /**
307 * Lab frame momentum (Px, Py, Pz, E and M in GeV) for the particle
308 * entries in this event.
309 */
310 std::vector< std::vector<double> > PUP;
311
312 /**
313 * Invariant lifetime (c*tau, distance from production to decay in
314 * mm) for the particle entries in this event.
315 */
316 std::vector<double> VTIMUP;
317
318 /**
319 * Spin info for the particle entries in this event given as the
320 * cosine of the angle between the spin vector of a particle and the
321 * 3-momentum of the decaying particle, specified in the lab frame.
322 */
323 std::vector<double> SPINUP;
324
325};
326
327/**
328 * The Reader class is initialized with a stream from which to read a
329 * version 1.0 Les Houches Accord event file. In the constructor of
330 * the Reader object the optional header information is read and then
331 * the mandatory init is read. After this the whole header block
332 * including the enclosing lines with tags are available in the public
333 * headerBlock member variable. Also the information from the init
334 * block is available in the heprup member variable and any additional
335 * comment lines are available in initComments. After each successful
336 * call to the readEvent() function the standard Les Houches Accord
337 * information about the event is available in the hepeup member
338 * variable and any additional comments in the eventComments
339 * variable. A typical reading sequence would look as follows:
340 *
341 *
342 */
343class Reader {
344
345public:
346
347 /**
348 * Initialize the Reader with a stream from which to read an event
349 * file. After the constructor is called the whole header block
350 * including the enclosing lines with tags are available in the
351 * public headerBlock member variable. Also the information from the
352 * init block is available in the heprup member variable and any
353 * additional comment lines are available in initComments.
354 *
355 * @param is the stream to read from.
356 */
357 Reader(std::istream & is)
358 : file(is) {
359 init();
360 }
361
362 /**
363 * Initialize the Reader with a filename from which to read an event
364 * file. After the constructor is called the whole header block
365 * including the enclosing lines with tags are available in the
366 * public headerBlock member variable. Also the information from the
367 * init block is available in the heprup member variable and any
368 * additional comment lines are available in initComments.
369 *
370 * @param filename the name of the file to read from.
371 */
372 Reader(std::string filename)
373 : intstream(filename.c_str()), file(intstream) {
374 init();
375 }
376
377private:
378
379 /**
380 * Used internally in the constructors to read header and init
381 * blocks.
382 */
383 void init() {
384
385 bool readingHeader = false;
386 bool readingInit = false;
387
388 // Make sure we are reading a LHEF file:
389 getline();
390 if ( currentLine.find("<LesHouchesEvents" ) == std::string::npos )
391 throw std::runtime_error
392 ("Tried to read a file which does not start with the "
393 "LesHouchesEvents tag.");
394 if ( currentLine.find("version=\"1.0\"" ) == std::string::npos )
395 throw std::runtime_error
396 ("Tried to read a LesHouchesEvents file which is not version 1.0.");
397
398 // Loop over all lines until we hit the </init> tag.
399 while ( getline() && currentLine.find("</init>") == std::string::npos ) {
400 if ( currentLine.find("<header") != std::string::npos ) {
401 // We have hit the header block, so we should dump this all
402 // following lines to headerBlock until we hit the end of it.
403 readingHeader = true;
404 headerBlock = currentLine + "\n";
405 }
406 else if ( currentLine.find("<init") != std::string::npos ) {
407 // We have hit the init block, so we should expect to find the
408 // standard information in the following.
409 readingInit = true;
410
411 // The first line tells us how many lines to read next.
412 getline();
413 std::istringstream iss(currentLine);
414 if ( !( iss >> heprup.IDBMUP.first >> heprup.IDBMUP.second
415 >> heprup.EBMUP.first >> heprup.EBMUP.second
416 >> heprup.PDFGUP.first >> heprup.PDFGUP.second
417 >> heprup.PDFSUP.first >> heprup.PDFSUP.second
418 >> heprup.IDWTUP >> heprup.NPRUP ) ) {
419 heprup.NPRUP = -42;
420 return;
421 }
422 heprup.resize();
423
424 for ( int i = 0; i < heprup.NPRUP; ++i ) {
425 getline();
426 std::istringstream iss(currentLine);
427 if ( !( iss >> heprup.XSECUP[i] >> heprup.XERRUP[i]
428 >> heprup.XMAXUP[i] >> heprup.LPRUP[i] ) ) {
429 heprup.NPRUP = -42;
430 return;
431 }
432 }
433 }
434 else if ( currentLine.find("</header>") != std::string::npos ) {
435 // The end of the header block. Dump this line as well to the
436 // headerBlock and we're done.
437 readingHeader = false;
438 headerBlock += currentLine + "\n";
439 }
440 else if ( readingHeader ) {
441 // We are in the process of reading the header block. Dump the
442 // line to haderBlock.
443 headerBlock += currentLine + "\n";
444 }
445 else if ( readingInit ) {
446 // Here we found a comment line. Dump it to initComments.
447 initComments += currentLine + "\n";
448 }
449 else {
450 // We found some other stuff outside the standard tags.
451 outsideBlock += currentLine + "\n";
452 }
453 }
454 if ( !file ) heprup.NPRUP = -42;
455 }
456
457public:
458
459 int getNumberOfEvents()
460 {
461 int counter = 0;
462 int position = file.tellg();
463 file.seekg(0, std::ios::beg);
464
465 file.clear();
466
467 while(getline())
468 {
469 if(currentLine.find("<event") != std::string::npos) ++counter;
470 }
471
472 file.clear();
473
474 file.seekg(position, std::ios::beg);
475
476 return counter;
477 }
478
479 /**
480 * Read an event from the file and store it in the hepeup
481 * object. Optional comment lines are stored i the eventComments
482 * member variable.
483 * @return true if the read sas successful.
484 */
485 bool readEvent() {
486
487 // Check if the initialization was successful. Otherwise we will
488 // not read any events.
489 if ( heprup.NPRUP < 0 ) return false;
490 eventComments = "";
491 outsideBlock = "";
492 hepeup.NUP = 0;
493
494 // Keep reading lines until we hit the next event or the end of
495 // the event block. Save any inbetween lines. Exit if we didn't
496 // find an event.
497 while ( getline() && currentLine.find("<event") == std::string::npos )
498 outsideBlock += currentLine + "\n";
499 if ( !getline() ) return false;
500
501 // We found an event. The first line determines how many
502 // subsequent particle lines we have.
503 std::istringstream iss(currentLine);
504 if ( !( iss >> hepeup.NUP >> hepeup.IDPRUP >> hepeup.XWGTUP
505 >> hepeup.SCALUP >> hepeup.AQEDUP >> hepeup.AQCDUP ) )
506 return false;
507 hepeup.resize();
508
509 // Read all particle lines.
510 for ( int i = 0; i < hepeup.NUP; ++i ) {
511 if ( !getline() ) return false;
512 std::istringstream iss(currentLine);
513 if ( !( iss >> hepeup.IDUP[i] >> hepeup.ISTUP[i]
514 >> hepeup.MOTHUP[i].first >> hepeup.MOTHUP[i].second
515 >> hepeup.ICOLUP[i].first >> hepeup.ICOLUP[i].second
516 >> hepeup.PUP[i][0] >> hepeup.PUP[i][1] >> hepeup.PUP[i][2]
517 >> hepeup.PUP[i][3] >> hepeup.PUP[i][4]
518 >> hepeup.VTIMUP[i] >> hepeup.SPINUP[i] ) )
519 return false;
520 }
521
522 // Now read any additional comments.
523 while ( getline() && currentLine.find("</event>") == std::string::npos )
524 eventComments += currentLine + "\n";
525
526 if ( !file ) return false;
527 return true;
528
529 }
530
531protected:
532
533 /**
534 * Used internally to read a single line from the stream.
535 */
536 bool getline() {
537 return ( std::getline(file, currentLine) );
538 }
539
540protected:
541
542 /**
543 * A local stream which is unused if a stream is supplied from the
544 * outside.
545 */
546 std::ifstream intstream;
547
548 /**
549 * The stream we are reading from. This may be a reference to an
550 * external stream or the internal intstream.
551 */
552 std::istream & file;
553
554 /**
555 * The last line read in from the stream in getline().
556 */
557 std::string currentLine;
558
559public:
560
561 /**
562 * All lines (since the last readEvent()) outside the header, init
563 * and event tags.
564 */
565 std::string outsideBlock;
566
567 /**
568 * All lines from the header block.
569 */
570 std::string headerBlock;
571
572 /**
573 * The standard init information.
574 */
575 HEPRUP heprup;
576
577 /**
578 * Additional comments found in the init block.
579 */
580 std::string initComments;
581
582 /**
583 * The standard information about the last read event.
584 */
585 HEPEUP hepeup;
586
587 /**
588 * Additional comments found with the last read event.
589 */
590 std::string eventComments;
591
592private:
593
594 /**
595 * The default constructor should never be used.
596 */
597 Reader();
598
599 /**
600 * The copy constructor should never be used.
601 */
602 Reader(const Reader &);
603
604 /**
605 * The Reader cannot be assigned to.
606 */
607 Reader & operator=(const Reader &);
608
609};
610
611/**
612 * The Writer class is initialized with a stream to which to write a
613 * version 1.0 Les Houches Accord event file. In the constructor of
614 * the Writer object the main XML tag is written out, with the
615 * corresponding end tag is written in the destructor. After a Writer
616 * object has been created, it is possible to assign standard init
617 * information in the heprup member variable. In addition any XML
618 * formatted information can be added to the headerBlock member
619 * variable (directly or via the addHeader() function). Further
620 * comment line (beginning with a <code>#</code> character) can be
621 * added to the initComments variable (directly or with the
622 * addInitComment() function). After this information is set, it
623 * should be written out to the file with the init() function.
624 *
625 * Before each event is written out with the writeEvent() function,
626 * the standard event information can then be assigned to the hepeup
627 * variable and optional comment lines (beginning with a
628 * <code>#</code> character) may be given to the eventComments
629 * variable (directly or with the addEventComment() function).
630 *
631 */
632class Writer {
633
634public:
635
636 /**
637 * Create a Writer object giving a stream to write to.
638 * @param os the stream where the event file is written.
639 */
640 Writer(std::ostream & os)
641 : file(os) {
642 // Write out the standard XML tag for the event file.
643 file << "<LesHouchesEvents version=\"1.0\">\n";
644 }
645
646 /**
647 * Create a Writer object giving a filename to write to.
648 * @param filename the name of the event file to be written.
649 */
650 Writer(std::string filename)
651 : intstream(filename.c_str()), file(intstream) {
652 // Write out the standard XML tag for the event file.
653 file << "LesHouchesEvents version=\"1.0\">\n";
654 }
655
656 /**
657 * The destructor writes out the final XML end-tag.
658 */
659 ~Writer() {
660 file << "</LesHouchesEvents>" << std::endl;
661 }
662
663 /**
664 * Add header lines consisting of XML code with this stream.
665 */
666 std::ostream & headerBlock() {
667 return headerStream;
668 }
669
670 /**
671 * Add comment lines to the init block with this stream.
672 */
673 std::ostream & initComments() {
674 return initStream;
675 }
676
677 /**
678 * Add comment lines to the next event to be written out with this stream.
679 */
680 std::ostream & eventComments() {
681 return eventStream;
682 }
683
684 /**
685 * Write out an optional header block followed by the standard init
686 * block information together with any comment lines.
687 */
688 void init() {
689
690 file << std::setprecision(8);
691
692 using std::setw;
693
694 std::string headerBlock = headerStream.str();
695 if ( headerBlock.length() ) {
696 if ( headerBlock.find("<header>") == std::string::npos )
697 file << "<header>\n";
698 if ( headerBlock[headerBlock.length() - 1] != '\n' )
699 headerBlock += '\n';
700 file << headerBlock;
701 if ( headerBlock.find("</header>") == std::string::npos )
702 file << "</header>\n";
703 }
704 file << "<init>\n"
705 << " " << setw(8) << heprup.IDBMUP.first
706 << " " << setw(8) << heprup.IDBMUP.second
707 << " " << setw(14) << heprup.EBMUP.first
708 << " " << setw(14) << heprup.EBMUP.second
709 << " " << setw(4) << heprup.PDFGUP.first
710 << " " << setw(4) << heprup.PDFGUP.second
711 << " " << setw(4) << heprup.PDFSUP.first
712 << " " << setw(4) << heprup.PDFSUP.second
713 << " " << setw(4) << heprup.IDWTUP
714 << " " << setw(4) << heprup.NPRUP << std::endl;
715 heprup.resize();
716 for ( int i = 0; i < heprup.NPRUP; ++i )
717 file << " " << setw(14) << heprup.XSECUP[i]
718 << " " << setw(14) << heprup.XERRUP[i]
719 << " " << setw(14) << heprup.XMAXUP[i]
720 << " " << setw(6) << heprup.LPRUP[i] << std::endl;
721 file << hashline(initStream.str()) << "</init>" << std::endl;
722 eventStream.str("");
723 }
724
725 /**
726 * Write out the event stored in hepeup, followed by optional
727 * comment lines.
728 */
729 bool writeEvent() {
730
731 using std::setw;
732
733 file << "<event>\n";
734 file << " " << setw(4) << hepeup.NUP
735 << " " << setw(6) << hepeup.IDPRUP
736 << " " << setw(14) << hepeup.XWGTUP
737 << " " << setw(14) << hepeup.SCALUP
738 << " " << setw(14) << hepeup.AQEDUP
739 << " " << setw(14) << hepeup.AQCDUP << "\n";
740 hepeup.resize();
741
742 for ( int i = 0; i < hepeup.NUP; ++i )
743 file << " " << setw(8) << hepeup.IDUP[i]
744 << " " << setw(2) << hepeup.ISTUP[i]
745 << " " << setw(4) << hepeup.MOTHUP[i].first
746 << " " << setw(4) << hepeup.MOTHUP[i].second
747 << " " << setw(4) << hepeup.ICOLUP[i].first
748 << " " << setw(4) << hepeup.ICOLUP[i].second
749 << " " << setw(14) << hepeup.PUP[i][0]
750 << " " << setw(14) << hepeup.PUP[i][1]
751 << " " << setw(14) << hepeup.PUP[i][2]
752 << " " << setw(14) << hepeup.PUP[i][3]
753 << " " << setw(14) << hepeup.PUP[i][4]
754 << " " << setw(1) << hepeup.VTIMUP[i]
755 << " " << setw(1) << hepeup.SPINUP[i] << std::endl;
756
757 file << hashline(eventStream.str()) << "</event>\n";
758
759 eventStream.str("");
760
761 if ( !file ) return false;
762
763 return true;
764
765 }
766
767protected:
768
769 /**
770 * Make sure that each line in the string \a s starts with a
771 * #-character and that the string ends with a new-line.
772 */
773 std::string hashline(std::string s) {
774 std::string ret;
775 std::istringstream is(s);
776 std::string ss;
777 while ( getline(is, ss) ) {
778 if ( ss.find('#') == std::string::npos ||
779 ss.find('#') != ss.find_first_not_of(" \t") ) ss = "# " + ss;
780 ret += ss + '\n';
781 }
782 return ret;
783 }
784
785protected:
786
787 /**
788 * A local stream which is unused if a stream is supplied from the
789 * outside.
790 */
791 std::ofstream intstream;
792
793 /**
794 * The stream we are writing to. This may be a reference to an
795 * external stream or the internal intstream.
796 */
797 std::ostream & file;
798
799public:
800
801 /**
802 * Stream to add all lines in the header block.
803 */
804 std::ostringstream headerStream;
805
806 /**
807 * The standard init information.
808 */
809 HEPRUP heprup;
810
811 /**
812 * Stream to add additional comments to be put in the init block.
813 */
814 std::ostringstream initStream;
815
816 /**
817 * The standard information about the event we will write next.
818 */
819 HEPEUP hepeup;
820
821 /**
822 * Stream to add additional comments to be written together the next event.
823 */
824 std::ostringstream eventStream;
825
826private:
827
828 /**
829 * The default constructor should never be used.
830 */
831 Writer();
832
833 /**
834 * The copy constructor should never be used.
835 */
836 Writer(const Writer &);
837
838 /**
839 * The Writer cannot be assigned to.
840 */
841 Writer & operator=(const Writer &);
842
843};
844
845}
846
847/** \example LHEFReadEx.cc An example function which reads from a Les
848 Huches Event File: */
849/** \example LHEFWriteEx.cc An example function which writes out a Les
850 Huches Event File: */
851/** \example LHEFCat.cc This is a main function which simply reads a
852 Les Houches Event File from the standard input and writes it again
853 to the standard output.
854 This file can be downloaded from
855 <A HREF="http://www.thep.lu.se/~leif/LHEF/LHEFCat.cc">here</A>.
856 There is also a sample
857 <A HREF="http://www.thep.lu.se/~leif/LHEF/ttbar.lhef">event file</A>
858 to try it on.
859*/
860
861/**\mainpage Les Houches Event File
862
863Why isn't any doxygen output generated by this text?
864
865Here are some example classes for reading and writing Les Houches
866Event Files according to the
867<A HREF="http://www.thep.lu.se/~torbjorn/lhef/lhafile2.pdf">proposal</A>
868by Torbj&ouml;rn Sj&ouml;strand discussed at the
869<A HREF="http://mc4lhc06.web.cern.ch/mc4lhc06/">MC4LHC</A>
870workshop at CERN 2006.
871
872In total there are four classes which are all available in a single
873header file called
874<A HREF="http://www.thep.lu.se/~leif/LHEF/LHEF.h">LHEF.h</A>.
875
876The two classes LHEF::HEPRUP and LHEF::HEPEUP are simple container
877classes which contain the same information as the Les Houches standard
878Fortran common blocks with the same names. The other two classes are
879called LHEF::Reader and LHEF::Writer and are used to read and write
880Les Houches Event Files
881
882Here are a few <A HREF="examples.html">examples</A> of how to use the
883classes:
884
885\namespace LHEF The LHEF namespace contains some example classes for reading and writing Les Houches
886Event Files according to the
887<A HREF="http://www.thep.lu.se/~torbjorn/lhef/lhafile2.pdf">proposal</A>
888by Torbj&ouml;rn Sj&ouml;strand discussed at the
889<A HREF="http://mc4lhc06.web.cern.ch/mc4lhc06/">MC4LHC</A>
890workshop at CERN 2006.
891
892
893
894 */
895
896
897#endif /* THEPEG_LHEF_H */
Note: See TracBrowser for help on using the repository browser.