Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooGenFitStudy.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 /**
18 \file RooGenFitStudy.cxx
19 \class RooGenFitStudy
20 \ingroup Roofitcore
21 
22 RooGenFitStudy is an abstract base class for RooStudyManager modules
23 
24 **/
25 
26 
27 
28 #include "RooFit.h"
29 #include "Riostream.h"
30 
31 #include "RooGenFitStudy.h"
32 #include "RooWorkspace.h"
33 #include "RooMsgService.h"
34 #include "RooDataSet.h"
35 #include "RooAbsPdf.h"
36 #include "RooRealVar.h"
37 #include "RooGlobalFunc.h"
38 #include "RooFitResult.h"
39 
40 
41 using namespace std ;
42 
43 ClassImp(RooGenFitStudy);
44  ;
45 
46 
47 ////////////////////////////////////////////////////////////////////////////////
48 /// Constructor
49 
50 RooGenFitStudy::RooGenFitStudy(const char* name, const char* title) :
51  RooAbsStudy(name?name:"RooGenFitStudy",title?title:"RooGenFitStudy"),
52  _genPdf(0),
53  _fitPdf(0),
54  _genSpec(0),
55  _nllVar(0),
56  _ngenVar(0),
57  _params(0),
58  _initParams(0)
59 {
60 }
61 
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 /// Copy constructor
66 
67 RooGenFitStudy::RooGenFitStudy(const RooGenFitStudy& other) :
68  RooAbsStudy(other),
69  _genPdfName(other._genPdfName),
70  _genObsName(other._genObsName),
71  _fitPdfName(other._fitPdfName),
72  _fitObsName(other._fitObsName),
73  _genPdf(0),
74  _fitPdf(0),
75  _genSpec(0),
76  _nllVar(0),
77  _ngenVar(0),
78  _params(0),
79  _initParams(0)
80 {
81  TIterator* giter = other._genOpts.MakeIterator() ;
82  TObject* o ;
83  while((o=giter->Next())) {
84  _genOpts.Add(o->Clone()) ;
85  }
86  delete giter ;
87 
88  TIterator* fiter = other._fitOpts.MakeIterator() ;
89  while((o=fiter->Next())) {
90  _fitOpts.Add(o->Clone()) ;
91  }
92  delete fiter ;
93 
94 }
95 
96 
97 
98 ////////////////////////////////////////////////////////////////////////////////
99 
100 RooGenFitStudy::~RooGenFitStudy()
101 {
102  if (_params) delete _params ;
103 }
104 
105 
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 /// Function called after insertion into workspace
109 
110 Bool_t RooGenFitStudy::attach(RooWorkspace& w)
111 {
112  Bool_t ret = kFALSE ;
113 
114  RooAbsPdf* pdf = w.pdf(_genPdfName.c_str()) ;
115  if (pdf) {
116  _genPdf = pdf ;
117  } else {
118  coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: generator p.d.f named " << _genPdfName << " not found in workspace " << w.GetName() << endl ;
119  ret = kTRUE ;
120  }
121 
122  _genObs.add(w.argSet(_genObsName.c_str())) ;
123  if (_genObs.getSize()==0) {
124  coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: no generator observables defined" << endl ;
125  ret = kTRUE ;
126  }
127 
128  pdf = w.pdf(_fitPdfName.c_str()) ;
129  if (pdf) {
130  _fitPdf = pdf ;
131  } else {
132  coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: fitting p.d.f named " << _fitPdfName << " not found in workspace " << w.GetName() << endl ;
133  ret = kTRUE ;
134  }
135 
136  _fitObs.add(w.argSet(_fitObsName.c_str())) ;
137  if (_fitObs.getSize()==0) {
138  coutE(InputArguments) << "RooGenFitStudy(" << GetName() << ") ERROR: no fitting observables defined" << endl ;
139  ret = kTRUE ;
140  }
141 
142  return ret ;
143 }
144 
145 
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 
149 void RooGenFitStudy::setGenConfig(const char* pdfName, const char* obsName, const RooCmdArg& arg1,const RooCmdArg& arg2,const RooCmdArg& arg3)
150 {
151  _genPdfName = pdfName ;
152  _genObsName = obsName ;
153  _genOpts.Add(arg1.Clone()) ;
154  _genOpts.Add(arg2.Clone()) ;
155  _genOpts.Add(arg3.Clone()) ;
156 }
157 
158 
159 
160 ////////////////////////////////////////////////////////////////////////////////
161 
162 void RooGenFitStudy::setFitConfig(const char* pdfName, const char* obsName, const RooCmdArg& arg1,const RooCmdArg& arg2,const RooCmdArg& arg3)
163 {
164  _fitPdfName = pdfName ;
165  _fitObsName = obsName ;
166  _fitOpts.Add(arg1.Clone()) ;
167  _fitOpts.Add(arg2.Clone()) ;
168  _fitOpts.Add(arg3.Clone()) ;
169 }
170 
171 
172 
173 ////////////////////////////////////////////////////////////////////////////////
174 /// One-time initialization of study
175 
176 Bool_t RooGenFitStudy::initialize()
177 {
178  _nllVar = new RooRealVar("NLL","-log(Likelihood)",0) ;
179  _ngenVar = new RooRealVar("ngen","number of generated events",0) ;
180 
181  _params = _fitPdf->getParameters(_genObs) ;
182  RooArgSet modelParams(*_params) ;
183  _initParams = (RooArgSet*) _params->snapshot() ;
184  _params->add(*_nllVar) ;
185  _params->add(*_ngenVar) ;
186 
187  _genSpec = _genPdf->prepareMultiGen(_genObs,(RooCmdArg&)*_genOpts.At(0),(RooCmdArg&)*_genOpts.At(1),(RooCmdArg&)*_genOpts.At(2)) ;
188 
189  registerSummaryOutput(*_params,modelParams) ;
190  return kFALSE ;
191 }
192 
193 
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 /// Execute one study iteration
197 
198 Bool_t RooGenFitStudy::execute()
199 {
200  *_params = *_initParams ;
201  RooDataSet* data = _genPdf->generate(*_genSpec) ;
202  RooFitResult* fr = _fitPdf->fitTo(*data,RooFit::Save(kTRUE),(RooCmdArg&)*_fitOpts.At(0),(RooCmdArg&)*_fitOpts.At(1),(RooCmdArg&)*_fitOpts.At(2)) ;
203 
204  if (fr->status()==0) {
205  _ngenVar->setVal(data->sumEntries()) ;
206  _nllVar->setVal(fr->minNll()) ;
207  storeSummaryOutput(*_params) ;
208  storeDetailedOutput(*fr) ;
209  }
210 
211  delete data ;
212  return kFALSE ;
213 }
214 
215 
216 
217 ////////////////////////////////////////////////////////////////////////////////
218 /// Finalization of study
219 
220 Bool_t RooGenFitStudy::finalize()
221 {
222  delete _params ;
223  delete _nllVar ;
224  delete _ngenVar ;
225  delete _initParams ;
226  delete _genSpec ;
227  _params = 0 ;
228  _nllVar = 0 ;
229  _ngenVar = 0 ;
230  _initParams = 0 ;
231  _genSpec = 0 ;
232 
233 
234  return kFALSE ;
235 }
236 
237 
238 ////////////////////////////////////////////////////////////////////////////////
239 
240 void RooGenFitStudy::Print(Option_t* /*options*/) const
241 {
242 }
243 
244