MINT2
Public Member Functions | Protected Attributes | Static Private Attributes | List of all members
DalitzEventList Class Reference

#include <DalitzEventList.h>

Inheritance diagram for DalitzEventList:
MINT::EventList< DalitzEvent > MINT::IEventList< DalitzEvent > MINT::PolymorphVector< DalitzEvent > MINT::IMinimalEventList< DalitzEvent >

Public Member Functions

 DalitzEventList ()
 
 DalitzEventList (const DalitzEventList &other)
 
 DalitzEventList (TNtupleD *)
 
virtual ~DalitzEventList ()
 
const std::string & className () const
 
int generatePhaseSpaceEvents (int NumEvents, const DalitzEventPattern &pat, TRandom *rnd=0)
 
int generateEvents (unsigned int NumEvents, const DalitzEventPattern &pat, MINT::IReturnRealForEvent< IDalitzEvent > *amps, TRandom *rnd)
 
TH1D * makePlot (const std::vector< int > sij, const std::string &name, MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, int nbins=100, double units=GeV *GeV, char opt='s')
 
TH2D * makePlot2D (const std::vector< int > sijIndicesX, const std::vector< int > sijIndicesY, const std::string &name, MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, int nbins=10, double units=GeV *GeV, char opt='s')
 
TNtupleD * makePlotNtp (const std::string &name_prefix="DalitzPlotNtp", MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, double units=GeV *GeV)
 
DalitzHistoSet histoSet () const
 
DalitzHistoSet weightedHistoSet () const
 
DalitzHistoSet reWeightedHistoSet (MINT::IReturnRealForEvent< IDalitzEvent > *w)
 
DalitzHistoSet weighedReWeightedHistoSet (MINT::IReturnRealForEvent< IDalitzEvent > *w)
 
bool makePlots (const std::string &filename) const
 
PlotSet makeAllPlots (const std::string &name_prefix, MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, int nbins1D=100, int nbins2D=10, double units=GeV *GeV)
 
bool save (const std::string &fname="DalitzEvents.root") const
 
bool fromFile (const std::string &fname="DalitzEvents.root")
 
TNtupleD * makeNtuple () const
 
TNtupleD * makeNtuple (const std::string &ntpName) const
 
bool saveAsNtuple (const std::string &fname="DalitzEvents.root") const
 
bool saveAsNtuple (const std::string &fname, const std::string &ntpName) const
 
bool fromNtuple (TTree *ntp)
 
bool fromNtuple (TTree *ntp, double num)
 
bool fromNtupleFile (const std::string &fname="DalitzEvents.root")
 
- Public Member Functions inherited from MINT::EventList< DalitzEvent >
 EventList ()
 
 EventList (const EventList< DalitzEvent > &other)
 
virtual ~EventList ()
 
EventList< DalitzEvent > & operator= (const EventList< DalitzEvent > &other)
 
virtual DalitzEventoperator[] (unsigned int i)
 
virtual const DalitzEventoperator[] (unsigned int i) const
 
virtual DalitzEvent getEvent (unsigned int i) const
 
virtual unsigned int size () const
 
virtual bool Add (const DalitzEvent &evt)
 
virtual bool Add (const EventList< DalitzEvent > &addList)
 
- Public Member Functions inherited from MINT::IEventList< DalitzEvent >
virtual ~IEventList ()
 
- Public Member Functions inherited from MINT::PolymorphVector< DalitzEvent >
 PolymorphVector ()
 
 PolymorphVector (unsigned int N)
 
 PolymorphVector (unsigned int N, const DalitzEvent &c)
 
 PolymorphVector (const PolymorphVector &other)
 
 PolymorphVector (const typename std::vector< DalitzEvent > &other)
 
virtual ~PolymorphVector ()
 
std::vector< DalitzEvent > & theVector ()
 
const std::vector< DalitzEvent > & theVector () const
 
DalitzEventoperator[] (unsigned int i)
 
const DalitzEventoperator[] (unsigned int i) const
 
DalitzEventat (unsigned int i)
 
const DalitzEventat (unsigned int i) const
 
std::vector< DalitzEvent >::iterator begin ()
 
std::vector< DalitzEvent >::const_iterator begin () const
 
std::vector< DalitzEvent >::iterator end ()
 
std::vector< DalitzEvent >::const_iterator end () const
 
std::vector< DalitzEvent >::iterator find (const DalitzEvent &c)
 
std::vector< DalitzEvent >::const_iterator find (const DalitzEvent &c) const
 
DalitzEventfront ()
 
const DalitzEventfront () const
 
DalitzEventback ()
 
const DalitzEventback () const
 
unsigned int size () const
 
bool empty () const
 
void push_back (const DalitzEvent &c)
 
void pop_back ()
 
void erase (typename std::vector< DalitzEvent >::iterator pos)
 
void erase (typename std::vector< DalitzEvent >::iterator first, typename std::vector< DalitzEvent >::iterator last)
 
PolymorphVector< DalitzEvent > & operator= (const PolymorphVector< DalitzEvent > &other)
 
void clear ()
 
void resize (unsigned int N)
 
void resize (unsigned int N, const DalitzEvent &c)
 
 operator const typename std::vector< DalitzEvent > & () const
 
 operator typename std::vector< DalitzEvent > & ()
 
bool operator== (const MINT::PolymorphVector< DalitzEvent > &v2) const
 
bool operator!= (const MINT::PolymorphVector< DalitzEvent > &v2) const
 
bool operator< (const MINT::PolymorphVector< DalitzEvent > &v2) const
 
bool operator> (const MINT::PolymorphVector< DalitzEvent > &v2) const
 

Protected Attributes

TFile * _file
 
- Protected Attributes inherited from MINT::PolymorphVector< DalitzEvent >
std::vector< DalitzEvent_vec
 

Static Private Attributes

static const std::string _className
 

Additional Inherited Members

- Protected Member Functions inherited from MINT::IEventList< DalitzEvent >
 IEventList ()
 
- Protected Member Functions inherited from MINT::IMinimalEventList< DalitzEvent >
 IMinimalEventList ()
 

Detailed Description

Definition at line 35 of file DalitzEventList.h.

Constructor & Destructor Documentation

◆ DalitzEventList() [1/3]

DalitzEventList::DalitzEventList ( )

Definition at line 27 of file MyDalitzEventList.h.

◆ DalitzEventList() [2/3]

DalitzEventList::DalitzEventList ( const DalitzEventList other)

◆ DalitzEventList() [3/3]

DalitzEventList::DalitzEventList ( TNtupleD *  ntp)

Definition at line 38 of file MyDalitzEventList.h.

40 {
41  bool success = fromNtuple(ntp);
42  if(! success){
43  cout << "ERROR in DalitzEventList::DalitzEventList(TNtupleD* ntp)"
44  << "\n > problem creating myself from ntuple"
45  << endl;
46  }
47 }
bool fromNtuple(TTree *ntp)

◆ ~DalitzEventList()

DalitzEventList::~DalitzEventList ( )
virtual

Definition at line 54 of file DalitzEventList.cpp.

54  {
55  //if(0 != _file) _file->Close();
56 }

Member Function Documentation

◆ className()

const std::string& DalitzEventList::className ( ) const
inline

Definition at line 49 of file DalitzEventList.h.

49 {return _className;}
static const std::string _className

◆ fromFile()

bool DalitzEventList::fromFile ( const std::string &  fname = "DalitzEvents.root")

Definition at line 320 of file MyDalitzEventList.h.

320  {
321  return fromNtupleFile(fname);
322 }
bool fromNtupleFile(const std::string &fname="DalitzEvents.root")

◆ fromNtuple() [1/2]

bool DalitzEventList::fromNtuple ( TTree *  ntp)

Definition at line 352 of file MyDalitzEventList.h.

352  {
353  bool dbThis=false;
354  if(dbThis) cout << "about to read ntuple with ptr " << ntp << endl;
355  if(0==ntp) return false;
356  if(ntp->GetEntries() <=0) return false;
357  if(dbThis) cout << " number of entries: " << ntp->GetEntries() << endl;
358  //if(dbThis) cout << " number of variables " << ntp->GetNvar() << endl;
359  bool success=true;
360  for(Long64_t i=0; i< ntp->GetEntries(); i++){
361  if(dbThis){
362  cout << "DalitzEventList::fromNtuple "
363  << " getting " << i << " th entry" << endl;
364  }
365  ntp->GetEntry(i);
366  if(dbThis) cout << " got it" << endl;
367  DalitzEvent evt;
368  // success &= evt.fromNtuple(ntp);
369  success &= evt.fromTree(ntp);
370  if(dbThis) cout << " made event" << endl;
371  if(! success){
372  cout << "ERROR in DalitzEventList::fromNtuple"
373  << ", call to DalitzEvent::fromNtuple returned false!"
374  << endl;
375  return false;
376  }
377  this->Add(evt);
378  if(dbThis) cout << " added event" << endl;
379  }
380  if(dbThis) cout << "DalitzEventList::fromNtuple worked!" << endl;
381  return success;
382 }
virtual bool Add(const DalitzEvent &evt)
Definition: EventList.h:63
bool fromTree(TTree *tree)

◆ fromNtuple() [2/2]

bool DalitzEventList::fromNtuple ( TTree *  ntp,
double  num 
)

Definition at line 384 of file MyDalitzEventList.h.

384  {
385  bool dbThis=false;
386  if(dbThis) cout << "about to read ntuple with ptr " << ntp << endl;
387  if(0==ntp) return false;
388  if(ntp->GetEntries() <=0) return false;
389  if(dbThis) cout << " number of entries: " << ntp->GetEntries() << endl;
390  //if(dbThis) cout << " number of variables " << ntp->GetNvar() << endl;
391  TRandom Rand(500);
392  bool success=true;
393  for(Long64_t i=0; i< ntp->GetEntries(); i++){
394  if(dbThis){
395  cout << "DalitzEventList::fromNtuple "
396  << " getting " << i << " th entry" << endl;
397  }
398  double rand = Rand.Rndm();
399  if (rand < num)
400  {
401  ntp->GetEntry(i);
402  if(dbThis) cout << " got it" << endl;
403  DalitzEvent evt;
404  // success &= evt.fromNtuple(ntp);
405  success &= evt.fromTree(ntp);
406  if(dbThis) cout << " made event" << endl;
407  if(! success){
408  cout << "ERROR in DalitzEventList::fromNtuple"
409  << ", call to DalitzEvent::fromNtuple returned false!"
410  << endl;
411  return false;
412  }
413  this->Add(evt);
414  if(dbThis) cout << " added event" << endl;
415  }
416  }
417  if(dbThis) cout << "DalitzEventList::fromNtuple worked!" << endl;
418  return success;
419 }
virtual bool Add(const DalitzEvent &evt)
Definition: EventList.h:63
bool fromTree(TTree *tree)

◆ fromNtupleFile()

bool DalitzEventList::fromNtupleFile ( const std::string &  fname = "DalitzEvents.root")

Definition at line 421 of file MyDalitzEventList.h.

421  {
422  TFile f(fname.c_str());
423  f.cd();
424  TTree* ntp = (TTree*) f.Get(className().c_str());
425  if(0 == ntp){
426  cout << "ERROR in DalitzEventList::fromNtupleFile"
427  << "\n > Can't get ntuple for filename = "
428  << "\n > " << fname << endl;
429  return false;
430  }
431  return fromNtuple(ntp);
432 }
bool fromNtuple(TTree *ntp)
const std::string & className() const

◆ generateEvents()

int DalitzEventList::generateEvents ( unsigned int  NumEvents,
const DalitzEventPattern pat,
MINT::IReturnRealForEvent< IDalitzEvent > *  amps,
TRandom *  rnd 
)

◆ generatePhaseSpaceEvents()

int DalitzEventList::generatePhaseSpaceEvents ( int  NumEvents,
const DalitzEventPattern pat,
TRandom *  rnd = 0 
)

Definition at line 49 of file MyDalitzEventList.h.

52  {
53 
54  cout << "DalitzEventList::generatePhaseSpaceEvents "
55  << " about to generate " << NumEvents
56  << " events to phase space" << endl;
57  time_t tstart = time(0);
58  int reportN = NumEvents/10; // 10 x at least;
59  if(reportN > 10000) reportN = 10000; // at least once every 10k event
60  for(int i=0; i< NumEvents; i++){
61  if( i < 5 || (!(i%reportN))){
62  cout << "DalitzEventList::generatePhaseSpaceEvents "
63  << " generating " << i+1 << "th event."
64  << " with " << pat << " and " << rnd << endl;
65  }
66  Add(DalitzEvent(pat, rnd));
67  }
68  double delT = difftime(time(0), tstart);
69  cout << " This took " << delT << " s";
70  if(delT > 0){
71  cout << " this is " << NumEvents/delT << " evts/s"
72  << " or " << (NumEvents/delT) * 60<< " evts/min";
73  }
74  cout << endl;
75 
76  return this->size();
77 }
virtual bool Add(const DalitzEvent &evt)
Definition: EventList.h:63
virtual unsigned int size() const
Definition: EventList.h:59

◆ histoSet()

DalitzHistoSet DalitzEventList::histoSet ( ) const

Definition at line 434 of file MyDalitzEventList.h.

434  {
435  DalitzHistoSet hs;
436  for(unsigned int i=0; i< this->size(); i++){
437  hs.addEvent(((*this)[i]));
438  }
439  return hs;
440 }
void addEvent(const IDalitzEvent &evt, double weight=1)
virtual unsigned int size() const
Definition: EventList.h:59

◆ makeAllPlots()

PlotSet DalitzEventList::makeAllPlots ( const std::string &  name_prefix,
MINT::IReturnRealForEvent< IDalitzEvent > *  weightFunction = 0,
int  nbins1D = 100,
int  nbins2D = 10,
double  units = GeV*GeV 
)

Definition at line 223 of file MyDalitzEventList.h.

228  {
229  PlotSet ps;
230  if(this->empty()) return ps;
231 
232  int nd = ((*this)[0]).eventPattern().size()-1;
233  if(nd < 3) return ps;
234 
235  AllPossibleSij sijList(nd);
236 
237  for(namedVMap::const_iterator it = sijList.begin();
238  it != sijList.end(); it++){
239  std::string name = name_prefix + "_s" + it->first;
240 
241  TH1D* newPlots = makePlot(it->second
242  , name
243  , weightFunction
244  , nbins1D
245  , units
246  , 's');
247  if(0 != newPlots) ps.push_back(newPlots);
248 
249  name = name_prefix + "_m" + it->first;
250  TH1D* newPlotm = makePlot(it->second
251  , name
252  , weightFunction
253  , nbins1D
254  , units
255  , 'm');
256  if(0 != newPlotm) ps.push_back(newPlotm);
257  namedVMap::const_iterator jt = it;
258  jt++;
259  for(; jt != sijList.end(); jt++){
260  std::string name2D = name_prefix + jt->first + "_vs_" + it->first;
261 
262  TH2D* newPlot2D = makePlot2D(it->second, jt->second, name2D
263  , weightFunction, nbins2D, units
264  , 's');
265 
266  if(0 != newPlot2D) ps.push_back(newPlot2D);
267  }
268  }
269  TNtupleD* newPlotNtp = makePlotNtp(name_prefix, weightFunction, units);
270  ps.push_back(newPlotNtp);
271 
272  return ps;
273 }
void push_back(const T &c)
TH2D * makePlot2D(const std::vector< int > sijIndicesX, const std::vector< int > sijIndicesY, const std::string &name, MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, int nbins=10, double units=GeV *GeV, char opt='s')
TNtupleD * makePlotNtp(const std::string &name_prefix="DalitzPlotNtp", MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, double units=GeV *GeV)
unsigned int size() const
TH1D * makePlot(const std::vector< int > sij, const std::string &name, MINT::IReturnRealForEvent< IDalitzEvent > *weightFunction=0, int nbins=100, double units=GeV *GeV, char opt='s')

◆ makeNtuple() [1/2]

TNtupleD * DalitzEventList::makeNtuple ( ) const

Definition at line 276 of file MyDalitzEventList.h.

276  {
277  return makeNtuple(className());
278 }
TNtupleD * makeNtuple() const
const std::string & className() const

◆ makeNtuple() [2/2]

TNtupleD * DalitzEventList::makeNtuple ( const std::string &  ntpName) const

Definition at line 286 of file MyDalitzEventList.h.

286  {
287 
288  if(this->empty()) return (TNtupleD*) 0;
289  std::string varNameString= ((this->begin()))->makeNtupleVarnames();
290  TNtupleD* ntp = new TNtupleD(className().c_str()
291  , ntpName.c_str()
292  , varNameString.c_str());
293 
294  ntp->SetDirectory(0);
295 
296  int arraySize = ((this->begin()))->ntupleVarArraySize();
297  Double_t *array = new Double_t[arraySize];
298 
299 
300  for(vector<DalitzEvent>::const_iterator it = this->begin();
301  it != this->end(); it++){
302 
303  bool success = (it)->fillNtupleVarArray(array, arraySize);
304  if(! success){
305  cout << "ERROR in DalitzEventList::makeNtuple"
306  << ", call to DalitzEvent::fillNtupleVarArray"
307  << " returned failure"
308  << endl;
309  }else{
310  ntp->Fill(array);
311  }
312  }
313  delete[] array;
314  return ntp;
315 }
std::vector< DalitzEvent >::iterator end()
std::vector< DalitzEvent >::iterator begin()
const std::string & className() const

◆ makePlot()

TH1D * DalitzEventList::makePlot ( const std::vector< int >  sij,
const std::string &  name,
MINT::IReturnRealForEvent< IDalitzEvent > *  weightFunction = 0,
int  nbins = 100,
double  units = GeV*GeV,
char  opt = 's' 
)

Definition at line 79 of file MyDalitzEventList.h.

85  {
86  if(this->empty()) return (TH1D*)0;
87 
88  double eps=0.1;
89  double mi = ((*this)[0]).sijMin(sijIndices)/units;
90  double ma = ((*this)[0]).sijMax(sijIndices)/units;
91  if(opt=='m'){
92  mi = sqrt(mi);
93  ma = sqrt(ma);
94  }
95  mi -= (ma-mi)*eps;
96  ma += (ma-mi)*eps;
97 
98  TH1D* histo = new TH1D(name.c_str(), name.c_str(), nbins, mi, ma);
99 
100  for(unsigned int i = 0; i < this->size(); i++){
101  double w=1;
102  DalitzEvent& evt((*this)[i]);
103  if(0 != weightFunction) w = weightFunction->RealVal(evt);
104  double x= evt.sij(sijIndices)/units;
105  if(opt == 'm') x = sqrt(x);
106  histo->Fill(x, w);
107  }
108  return histo;
109 }
virtual double RealVal(EVENT_TYPE &evt)=0
virtual unsigned int size() const
Definition: EventList.h:59

◆ makePlot2D()

TH2D * DalitzEventList::makePlot2D ( const std::vector< int >  sijIndicesX,
const std::vector< int >  sijIndicesY,
const std::string &  name,
MINT::IReturnRealForEvent< IDalitzEvent > *  weightFunction = 0,
int  nbins = 10,
double  units = GeV*GeV,
char  opt = 's' 
)

Definition at line 111 of file MyDalitzEventList.h.

118  {
119  if(this->empty()) return (TH2D*)0;
120 
121 
122  double eps = 0.125;
123  double miX = ((*this)[0]).sijMin(sijIndicesX)/units;
124  double maX = ((*this)[0]).sijMax(sijIndicesX)/units;
125  double miY = ((*this)[0]).sijMin(sijIndicesY)/units;
126  double maY = ((*this)[0]).sijMax(sijIndicesY)/units;
127  if(opt=='m'){
128  miX = sqrt(miX);
129  maX = sqrt(maX);
130  miY = sqrt(miY);
131  maY = sqrt(maY);
132  }
133  miX -= (maX - miX)*eps;
134  maX += (maX - miX)*eps;
135  miY -= (maY - miY)*eps;
136  maY += (maY - miY)*eps;
137 
138  TH2D* histo = new TH2D(name.c_str()
139  , name.c_str()
140  , nbins, miX, maX
141  , nbins, miY, maY);
142  for(unsigned int i = 0; i < this->size(); i++){
143  double w=1;
144  DalitzEvent& evt((*this)[i]);
145  if(0 != weightFunction) w = weightFunction->RealVal(evt);
146  double x= evt.sij(sijIndicesX)/units;
147  double y= evt.sij(sijIndicesY)/units;
148  if(opt == 'm'){
149  x = sqrt(x);
150  y = sqrt(y);
151  }
152  histo->Fill(x, y, w);
153  }
154  return histo;
155 }
virtual double RealVal(EVENT_TYPE &evt)=0
virtual unsigned int size() const
Definition: EventList.h:59

◆ makePlotNtp()

TNtupleD * DalitzEventList::makePlotNtp ( const std::string &  name_prefix = "DalitzPlotNtp",
MINT::IReturnRealForEvent< IDalitzEvent > *  weightFunction = 0,
double  units = GeV*GeV 
)

Definition at line 157 of file MyDalitzEventList.h.

160  {
161  std::string ntp_str="";
162  int nd = ((*this)[0]).eventPattern().size()-1;
163  if(nd < 3) return (TNtupleD*)0;
164  AllPossibleSij sijList(nd);
165 
166  int arraySize = 0;
167  for(namedVMap::const_iterator it = sijList.begin();
168  it != sijList.end(); it++){
169  std::string entry_name = "s" + it->first + ":";
170  ntp_str += entry_name;
171  arraySize++;
172  }
173  for(namedVMap::const_iterator it = sijList.begin();
174  it != sijList.end(); it++){
175  std::string entry_name = "m" + it->first + ":";
176  ntp_str += entry_name;
177  arraySize++;
178  }
179  ntp_str += "weight:";
180  arraySize++;
181  ntp_str += "fcn"; // last one doesn't get a ":
182 
183  std::string ntp_id = className() + "_sij";
184  std::string name = name_prefix + "_sij";
185  TNtupleD* ntp = new TNtupleD(ntp_id.c_str()
186  , name.c_str()
187  , ntp_str.c_str());
188 
189  Double_t* array = new Double_t[arraySize];
190 
191  for(unsigned int i = 0; i < this->size(); i++){
192  int arrayIndex = 0;
193  for(namedVMap::const_iterator it = sijList.begin();
194  it != sijList.end(); it++){
195  array[arrayIndex]=((*this)[i]).sij(it->second)/units;
196  arrayIndex++;
197  }
198  for(namedVMap::const_iterator it = sijList.begin();
199  it != sijList.end(); it++){
200  std::string entry_name = "m" + it->first;
201  double sijtemp = ((*this)[i]).sij(it->second)/units;
202  if(sijtemp >=0){
203  array[arrayIndex]=sqrt(sijtemp);
204  }else{
205  array[arrayIndex] = -9999.;
206  }
207  arrayIndex++;
208  }
209  DalitzEvent& evt((*this)[i]);
210  array[arrayIndex] = evt.getWeight();
211  arrayIndex++;
212  double fcn=1;
213  if(0 != weightFunction) fcn = weightFunction->RealVal(evt);
214  array[arrayIndex] = fcn;
215  arrayIndex++;
216 
217  ntp->Fill(array);
218  }
219  delete[] array;
220 
221  return ntp;
222 }
virtual double RealVal(EVENT_TYPE &evt)=0
virtual unsigned int size() const
Definition: EventList.h:59
const std::string & className() const

◆ makePlots()

bool DalitzEventList::makePlots ( const std::string &  filename) const

Definition at line 471 of file MyDalitzEventList.h.

471  {
472  histoSet().save(filename);
473  return true;
474 }
DalitzHistoSet histoSet() const
bool save(const std::string &filename="DalitzHistos.root") const

◆ reWeightedHistoSet()

DalitzHistoSet DalitzEventList::reWeightedHistoSet ( MINT::IReturnRealForEvent< IDalitzEvent > *  w)

Definition at line 450 of file MyDalitzEventList.h.

450  {
451  // mainly for diagnostics
452  DalitzHistoSet hs;
453  if(0 == w) return hs;
454  for(unsigned int i=0; i< this->size(); i++){
455  DalitzEvent evt((*this)[i]);
456  hs.addEvent(evt, w->RealVal(evt));
457  }
458  return hs;
459 }
virtual double RealVal(EVENT_TYPE &evt)=0
void addEvent(const IDalitzEvent &evt, double weight=1)
virtual unsigned int size() const
Definition: EventList.h:59

◆ save()

bool DalitzEventList::save ( const std::string &  fname = "DalitzEvents.root") const

Definition at line 317 of file MyDalitzEventList.h.

317  {
318  return saveAsNtuple(fname);
319 }
bool saveAsNtuple(const std::string &fname="DalitzEvents.root") const

◆ saveAsNtuple() [1/2]

bool DalitzEventList::saveAsNtuple ( const std::string &  fname = "DalitzEvents.root") const

Definition at line 324 of file MyDalitzEventList.h.

325  {
326  return saveAsNtuple(fname, className());
327 }
const std::string & className() const
bool saveAsNtuple(const std::string &fname="DalitzEvents.root") const

◆ saveAsNtuple() [2/2]

bool DalitzEventList::saveAsNtuple ( const std::string &  fname,
const std::string &  ntpName 
) const

Definition at line 329 of file MyDalitzEventList.h.

331  {
332  if(this->empty()){
333  cout << "WARNING in DalitzEventList::saveAsNtuple!"
334  << "\n\tyou are trying to save me to the file: "
335  << fname
336  << "\n\tbut I have only " << this->size()
337  << " events."
338  << " I won't create the file."
339  << endl;
340  return false;
341  }
342  TFile f(fname.c_str(), "RECREATE");
343  f.cd();
344  TNtupleD* ntp = makeNtuple(ntpName);
345  ntp->Write();
346  f.Close();
347  //ntp->Delete("all");
348  //delete ntp;
349  return true;
350 }
TNtupleD * makeNtuple() const
virtual unsigned int size() const
Definition: EventList.h:59

◆ weighedReWeightedHistoSet()

DalitzHistoSet DalitzEventList::weighedReWeightedHistoSet ( MINT::IReturnRealForEvent< IDalitzEvent > *  w)

Definition at line 461 of file MyDalitzEventList.h.

461  {
462  DalitzHistoSet hs;
463  if(0 == w) return hs;
464  for(unsigned int i=0; i< this->size(); i++){
465  DalitzEvent evt((*this)[i]);
466  hs.addEvent(evt, w->RealVal(evt) * evt.getWeight());
467  }
468  return hs;
469 }
virtual double RealVal(EVENT_TYPE &evt)=0
void addEvent(const IDalitzEvent &evt, double weight=1)
virtual unsigned int size() const
Definition: EventList.h:59

◆ weightedHistoSet()

DalitzHistoSet DalitzEventList::weightedHistoSet ( ) const

Definition at line 441 of file MyDalitzEventList.h.

441  {
442  // mainly for diagnostics
443  DalitzHistoSet hs;
444  for(unsigned int i=0; i< this->size(); i++){
445  DalitzEvent evt((*this)[i]);
446  hs.addEvent(evt, evt.getWeight());
447  }
448  return hs;
449 }
void addEvent(const IDalitzEvent &evt, double weight=1)
virtual unsigned int size() const
Definition: EventList.h:59

Member Data Documentation

◆ _className

const std::string DalitzEventList::_className
staticprivate

Definition at line 38 of file DalitzEventList.h.

◆ _file

TFile* DalitzEventList::_file
protected

Definition at line 40 of file DalitzEventList.h.


The documentation for this class was generated from the following files: