Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
RooAbsTestStatistic.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 RooAbsTestStatistic.cxx
19 \class RooAbsTestStatistic
20 \ingroup Roofitcore
21 
22 RooAbsTestStatistic is the abstract base class for all test
23 statistics. Test statistics that evaluate the PDF at each data
24 point should inherit from the RooAbsOptTestStatistic class which
25 implements several generic optimizations that can be done for such
26 quantities.
27 
28 This test statistic base class organizes calculation of test
29 statistic values for RooSimultaneous PDF as a combination of test
30 statistic values for the PDF components of the simultaneous PDF and
31 organizes multi-processor parallel calculation of test statistic
32 values. For the latter, the test statistic value is calculated in
33 partitions in parallel executing processes and a posteriori
34 combined in the main thread.
35 **/
36 
37 
38 #include "RooFit.h"
39 #include "Riostream.h"
40 
41 #include "RooAbsTestStatistic.h"
42 #include "RooAbsPdf.h"
43 #include "RooSimultaneous.h"
44 #include "RooAbsData.h"
45 #include "RooArgSet.h"
46 #include "RooRealVar.h"
47 #include "RooNLLVar.h"
48 #include "RooRealMPFE.h"
49 #include "RooErrorHandler.h"
50 #include "RooMsgService.h"
51 #include "TTimeStamp.h"
52 #include "RooProdPdf.h"
53 #include "RooRealSumPdf.h"
54 #include <string>
55 
56 using namespace std;
57 
58 ClassImp(RooAbsTestStatistic);
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Default constructor
62 
63 RooAbsTestStatistic::RooAbsTestStatistic() :
64  _func(0), _data(0), _projDeps(0), _splitRange(0), _simCount(0),
65  _verbose(kFALSE), _init(kFALSE), _gofOpMode(Slave), _nEvents(0), _setNum(0),
66  _numSets(0), _extSet(0), _nGof(0), _gofArray(0), _nCPU(1), _mpfeArray(0),
67  _mpinterl(RooFit::BulkPartition), _doOffset(kFALSE), _offset(0),
68  _offsetCarry(0), _evalCarry(0)
69 {
70 }
71 
72 
73 
74 ////////////////////////////////////////////////////////////////////////////////
75 /// Create a test statistic from the given function and the data.
76 /// \param[in] name Name of the test statistic
77 /// \param[in] title Title (for plotting)
78 /// \param[in] real Function to be used for tests
79 /// \param[in] data Data to fit function to
80 /// \param[in] projDeps A set of projected observables
81 /// \param[in] rangeName Fit data only in range with given name
82 /// \param[in] addCoefRangeName If not null, all RooAddPdf components of `real` will be instructed to fix their fraction definitions to the given named range.
83 /// \param[in] nCPU If larger than one, the test statistic calculation will be parallelized over multiple processes.
84 /// By default the data is split with 'bulk' partitioning (each process calculates a contigious block of fraction 1/nCPU
85 /// of the data). For binned data this approach may be suboptimal as the number of bins with >0 entries
86 /// in each processing block many vary greatly thereby distributing the workload rather unevenly.
87 /// \param[in] interleave is set to true, the interleave partitioning strategy is used where each partition
88 /// i takes all bins for which (ibin % ncpu == i) which is more likely to result in an even workload.
89 /// \param[in] verbose Be more verbose.
90 /// \param[in] splitCutRange If true, a different rangeName constructed as rangeName_{catName} will be used
91 /// as range definition for each index state of a RooSimultaneous. This means that a different range can be defined
92 /// for each category such as
93 /// ```
94 /// myVariable.setRange("range_pi0", 135, 210);
95 /// myVariable.setRange("range_gamma", 50, 210);
96 /// ```
97 /// if the categories are called "pi0" and "gamma".
98 
99 RooAbsTestStatistic::RooAbsTestStatistic(const char *name, const char *title, RooAbsReal& real, RooAbsData& data,
100  const RooArgSet& projDeps, const char* rangeName, const char* addCoefRangeName,
101  Int_t nCPU, RooFit::MPSplit interleave, Bool_t verbose, Bool_t splitCutRange) :
102  RooAbsReal(name,title),
103  _paramSet("paramSet","Set of parameters",this),
104  _func(&real),
105  _data(&data),
106  _projDeps((RooArgSet*)projDeps.Clone()),
107  _rangeName(rangeName?rangeName:""),
108  _addCoefRangeName(addCoefRangeName?addCoefRangeName:""),
109  _splitRange(splitCutRange),
110  _simCount(1),
111  _verbose(verbose),
112  _nGof(0),
113  _gofArray(0),
114  _nCPU(nCPU),
115  _mpfeArray(0),
116  _mpinterl(interleave),
117  _doOffset(kFALSE),
118  _offset(0),
119  _offsetCarry(0),
120  _evalCarry(0)
121 {
122  // Register all parameters as servers
123  RooArgSet* params = real.getParameters(&data) ;
124  _paramSet.add(*params) ;
125  delete params ;
126 
127  if (_nCPU>1 || _nCPU==-1) {
128 
129  if (_nCPU==-1) {
130  _nCPU=1 ;
131  }
132 
133  _gofOpMode = MPMaster ;
134 
135  } else {
136 
137  // Determine if RooAbsReal is a RooSimultaneous
138  Bool_t simMode = dynamic_cast<RooSimultaneous*>(&real)?kTRUE:kFALSE ;
139 
140  if (simMode) {
141  _gofOpMode = SimMaster ;
142  } else {
143  _gofOpMode = Slave ;
144  }
145  }
146 
147  _setNum = 0 ;
148  _extSet = 0 ;
149  _numSets = 1 ;
150  _init = kFALSE ;
151  _nEvents = data.numEntries() ;
152 }
153 
154 
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// Copy constructor
158 
159 RooAbsTestStatistic::RooAbsTestStatistic(const RooAbsTestStatistic& other, const char* name) :
160  RooAbsReal(other,name),
161  _paramSet("paramSet","Set of parameters",this),
162  _func(other._func),
163  _data(other._data),
164  _projDeps((RooArgSet*)other._projDeps->Clone()),
165  _rangeName(other._rangeName),
166  _addCoefRangeName(other._addCoefRangeName),
167  _splitRange(other._splitRange),
168  _simCount(1),
169  _verbose(other._verbose),
170  _nGof(0),
171  _gofArray(0),
172  _gofSplitMode(other._gofSplitMode),
173  _nCPU(other._nCPU),
174  _mpfeArray(0),
175  _mpinterl(other._mpinterl),
176  _doOffset(other._doOffset),
177  _offset(other._offset),
178  _offsetCarry(other._offsetCarry),
179  _evalCarry(other._evalCarry)
180 {
181  // Our parameters are those of original
182  _paramSet.add(other._paramSet) ;
183 
184  if (_nCPU>1 || _nCPU==-1) {
185 
186  if (_nCPU==-1) {
187  _nCPU=1 ;
188  }
189 
190  _gofOpMode = MPMaster ;
191 
192  } else {
193 
194  // Determine if RooAbsReal is a RooSimultaneous
195  Bool_t simMode = dynamic_cast<RooSimultaneous*>(_func)?kTRUE:kFALSE ;
196 
197  if (simMode) {
198  _gofOpMode = SimMaster ;
199  } else {
200  _gofOpMode = Slave ;
201  }
202  }
203 
204  _setNum = 0 ;
205  _extSet = 0 ;
206  _numSets = 1 ;
207  _init = kFALSE ;
208  _nEvents = _data->numEntries() ;
209 
210 
211 }
212 
213 
214 
215 ////////////////////////////////////////////////////////////////////////////////
216 /// Destructor
217 
218 RooAbsTestStatistic::~RooAbsTestStatistic()
219 {
220  if (MPMaster == _gofOpMode && _init) {
221  for (Int_t i = 0; i < _nCPU; ++i) delete _mpfeArray[i];
222  delete[] _mpfeArray ;
223  }
224 
225  if (SimMaster == _gofOpMode && _init) {
226  for (Int_t i = 0; i < _nGof; ++i) delete _gofArray[i];
227  delete[] _gofArray ;
228  }
229 
230  delete _projDeps ;
231 
232 }
233 
234 
235 
236 ////////////////////////////////////////////////////////////////////////////////
237 /// Calculate and return value of test statistic. If the test statistic
238 /// is calculated from a RooSimultaneous, the test statistic calculation
239 /// is performed separately on each simultaneous p.d.f component and associated
240 /// data, and then combined. If the test statistic calculation is parallelized,
241 /// partitions are calculated in nCPU processes and combined a posteriori.
242 
243 Double_t RooAbsTestStatistic::evaluate() const
244 {
245  // One-time Initialization
246  if (!_init) {
247  const_cast<RooAbsTestStatistic*>(this)->initialize() ;
248  }
249 
250  if (SimMaster == _gofOpMode) {
251  // Evaluate array of owned GOF objects
252  Double_t ret = 0.;
253 
254  if (_mpinterl == RooFit::BulkPartition || _mpinterl == RooFit::Interleave ) {
255  ret = combinedValue((RooAbsReal**)_gofArray,_nGof);
256  } else {
257  Double_t sum = 0., carry = 0.;
258  for (Int_t i = 0 ; i < _nGof; ++i) {
259  if (i % _numSets == _setNum || (_mpinterl==RooFit::Hybrid && _gofSplitMode[i] != RooFit::SimComponents )) {
260  Double_t y = _gofArray[i]->getValV();
261  carry += _gofArray[i]->getCarry();
262  y -= carry;
263  const Double_t t = sum + y;
264  carry = (t - sum) - y;
265  sum = t;
266  }
267  }
268  ret = sum ;
269  _evalCarry = carry;
270  }
271 
272  // Only apply global normalization if SimMaster doesn't have MP master
273  if (numSets()==1) {
274  const Double_t norm = globalNormalization();
275  ret /= norm;
276  _evalCarry /= norm;
277  }
278 
279  return ret ;
280 
281  } else if (MPMaster == _gofOpMode) {
282 
283  // Start calculations in parallel
284  for (Int_t i = 0; i < _nCPU; ++i) _mpfeArray[i]->calculate();
285 
286 
287  Double_t sum(0), carry = 0.;
288  for (Int_t i = 0; i < _nCPU; ++i) {
289  Double_t y = _mpfeArray[i]->getValV();
290  carry += _mpfeArray[i]->getCarry();
291  y -= carry;
292  const Double_t t = sum + y;
293  carry = (t - sum) - y;
294  sum = t;
295  }
296 
297  Double_t ret = sum ;
298  _evalCarry = carry;
299  return ret ;
300 
301  } else {
302 
303  // Evaluate as straight FUNC
304  Int_t nFirst(0), nLast(_nEvents), nStep(1) ;
305 
306  switch (_mpinterl) {
307  case RooFit::BulkPartition:
308  nFirst = _nEvents * _setNum / _numSets ;
309  nLast = _nEvents * (_setNum+1) / _numSets ;
310  nStep = 1 ;
311  break;
312 
313  case RooFit::Interleave:
314  nFirst = _setNum ;
315  nLast = _nEvents ;
316  nStep = _numSets ;
317  break ;
318 
319  case RooFit::SimComponents:
320  nFirst = 0 ;
321  nLast = _nEvents ;
322  nStep = 1 ;
323  break ;
324 
325  case RooFit::Hybrid:
326  throw(std::string("this should never happen")) ;
327  break ;
328  }
329 
330  Double_t ret = evaluatePartition(nFirst,nLast,nStep);
331 
332  if (numSets()==1) {
333  const Double_t norm = globalNormalization();
334  ret /= norm;
335  _evalCarry /= norm;
336  }
337 
338  return ret ;
339 
340  }
341 }
342 
343 
344 
345 ////////////////////////////////////////////////////////////////////////////////
346 /// One-time initialization of the test statistic. Setup
347 /// infrastructure for simultaneous p.d.f processing and/or
348 /// parallelized processing if requested
349 
350 Bool_t RooAbsTestStatistic::initialize()
351 {
352  if (_init) return kFALSE;
353 
354  if (MPMaster == _gofOpMode) {
355  initMPMode(_func,_data,_projDeps,_rangeName.size()?_rangeName.c_str():0,_addCoefRangeName.size()?_addCoefRangeName.c_str():0) ;
356  } else if (SimMaster == _gofOpMode) {
357  initSimMode((RooSimultaneous*)_func,_data,_projDeps,_rangeName.size()?_rangeName.c_str():0,_addCoefRangeName.size()?_addCoefRangeName.c_str():0) ;
358  }
359  _init = kTRUE;
360  return kFALSE;
361 }
362 
363 
364 
365 ////////////////////////////////////////////////////////////////////////////////
366 /// Forward server redirect calls to component test statistics
367 
368 Bool_t RooAbsTestStatistic::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t)
369 {
370  if (SimMaster == _gofOpMode && _gofArray) {
371  // Forward to slaves
372  for (Int_t i = 0; i < _nGof; ++i) {
373  if (_gofArray[i]) {
374  _gofArray[i]->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange);
375  }
376  }
377  } else if (MPMaster == _gofOpMode&& _mpfeArray) {
378  // Forward to slaves
379  for (Int_t i = 0; i < _nCPU; ++i) {
380  if (_mpfeArray[i]) {
381  _mpfeArray[i]->recursiveRedirectServers(newServerList,mustReplaceAll,nameChange);
382 // cout << "redirecting servers on " << _mpfeArray[i]->GetName() << endl;
383  }
384  }
385  }
386  return kFALSE;
387 }
388 
389 
390 
391 ////////////////////////////////////////////////////////////////////////////////
392 /// Add extra information on component test statistics when printing
393 /// itself as part of a tree structure
394 
395 void RooAbsTestStatistic::printCompactTreeHook(ostream& os, const char* indent)
396 {
397  if (SimMaster == _gofOpMode) {
398  // Forward to slaves
399  os << indent << "RooAbsTestStatistic begin GOF contents" << endl ;
400  for (Int_t i = 0; i < _nGof; ++i) {
401  if (_gofArray[i]) {
402  TString indent2(indent);
403  indent2 += Form("[%d] ",i);
404  _gofArray[i]->printCompactTreeHook(os,indent2);
405  }
406  }
407  os << indent << "RooAbsTestStatistic end GOF contents" << endl;
408  } else if (MPMaster == _gofOpMode) {
409  // WVE implement this
410  }
411 }
412 
413 
414 
415 ////////////////////////////////////////////////////////////////////////////////
416 /// Forward constant term optimization management calls to component
417 /// test statistics
418 
419 void RooAbsTestStatistic::constOptimizeTestStatistic(ConstOpCode opcode, Bool_t doAlsoTrackingOpt)
420 {
421  initialize();
422  if (SimMaster == _gofOpMode) {
423  // Forward to slaves
424  for (Int_t i = 0; i < _nGof; ++i) {
425  // In SimComponents Splitting strategy only constOptimize the terms that are actually used
426  RooFit::MPSplit effSplit = (_mpinterl!=RooFit::Hybrid) ? _mpinterl : _gofSplitMode[i];
427  if ( (i % _numSets == _setNum) || (effSplit != RooFit::SimComponents) ) {
428  if (_gofArray[i]) _gofArray[i]->constOptimizeTestStatistic(opcode,doAlsoTrackingOpt);
429  }
430  }
431  } else if (MPMaster == _gofOpMode) {
432  for (Int_t i = 0; i < _nCPU; ++i) {
433  _mpfeArray[i]->constOptimizeTestStatistic(opcode,doAlsoTrackingOpt);
434  }
435  }
436 }
437 
438 
439 
440 ////////////////////////////////////////////////////////////////////////////////
441 /// Set MultiProcessor set number identification of this instance
442 
443 void RooAbsTestStatistic::setMPSet(Int_t inSetNum, Int_t inNumSets)
444 {
445  _setNum = inSetNum; _numSets = inNumSets;
446  _extSet = _mpinterl==RooFit::SimComponents ? _setNum : (_numSets - 1);
447 
448  if (SimMaster == _gofOpMode) {
449  // Forward to slaves
450  initialize();
451  for (Int_t i = 0; i < _nGof; ++i) {
452  if (_gofArray[i]) _gofArray[i]->setMPSet(inSetNum,inNumSets);
453  }
454  }
455 }
456 
457 
458 
459 ////////////////////////////////////////////////////////////////////////////////
460 /// Initialize multi-processor calculation mode. Create component test statistics in separate
461 /// processed that are connected to this process through a RooAbsRealMPFE front-end class.
462 
463 void RooAbsTestStatistic::initMPMode(RooAbsReal* real, RooAbsData* data, const RooArgSet* projDeps, const char* rangeName, const char* addCoefRangeName)
464 {
465  _mpfeArray = new pRooRealMPFE[_nCPU];
466 
467  // Create proto-goodness-of-fit
468  RooAbsTestStatistic* gof = create(GetName(),GetTitle(),*real,*data,*projDeps,rangeName,addCoefRangeName,1,_mpinterl,_verbose,_splitRange);
469  gof->recursiveRedirectServers(_paramSet);
470 
471  for (Int_t i = 0; i < _nCPU; ++i) {
472  gof->setMPSet(i,_nCPU);
473  gof->SetName(Form("%s_GOF%d",GetName(),i));
474  gof->SetTitle(Form("%s_GOF%d",GetTitle(),i));
475 
476  ccoutD(Eval) << "RooAbsTestStatistic::initMPMode: starting remote server process #" << i << endl;
477  _mpfeArray[i] = new RooRealMPFE(Form("%s_%lx_MPFE%d",GetName(),(ULong_t)this,i),Form("%s_%lx_MPFE%d",GetTitle(),(ULong_t)this,i),*gof,false);
478  //_mpfeArray[i]->setVerbose(kTRUE,kTRUE);
479  _mpfeArray[i]->initialize();
480  if (i > 0) {
481  _mpfeArray[i]->followAsSlave(*_mpfeArray[0]);
482  }
483  }
484  _mpfeArray[_nCPU - 1]->addOwnedComponents(*gof);
485  coutI(Eval) << "RooAbsTestStatistic::initMPMode: started " << _nCPU << " remote server process." << endl;
486  //cout << "initMPMode --- done" << endl ;
487  return ;
488 }
489 
490 
491 
492 ////////////////////////////////////////////////////////////////////////////////
493 /// Initialize simultaneous p.d.f processing mode. Strip simultaneous
494 /// p.d.f into individual components, split dataset in subset
495 /// matching each component and create component test statistics for
496 /// each of them.
497 
498 void RooAbsTestStatistic::initSimMode(RooSimultaneous* simpdf, RooAbsData* data,
499  const RooArgSet* projDeps, const char* rangeName, const char* addCoefRangeName)
500 {
501 
502  RooAbsCategoryLValue& simCat = (RooAbsCategoryLValue&) simpdf->indexCat();
503 
504  TString simCatName(simCat.GetName());
505  TList* dsetList = const_cast<RooAbsData*>(data)->split(simCat,processEmptyDataSets());
506  if (!dsetList) {
507  coutE(Fitting) << "RooAbsTestStatistic::initSimMode(" << GetName() << ") ERROR: index category of simultaneous pdf is missing in dataset, aborting" << endl;
508  throw std::string("RooAbsTestStatistic::initSimMode() ERROR, index category of simultaneous pdf is missing in dataset, aborting");
509  //RooErrorHandler::softAbort() ;
510  }
511 
512  // Count number of used states
513  Int_t n = 0;
514  _nGof = 0;
515  RooCatType* type;
516  TIterator* catIter = simCat.typeIterator();
517  while ((type = (RooCatType*) catIter->Next())) {
518  // Retrieve the PDF for this simCat state
519  RooAbsPdf* pdf = simpdf->getPdf(type->GetName());
520  RooAbsData* dset = (RooAbsData*) dsetList->FindObject(type->GetName());
521 
522  if (pdf && dset && (0. != dset->sumEntries() || processEmptyDataSets())) {
523  ++_nGof;
524  }
525  }
526 
527  // Allocate arrays
528  _gofArray = new pRooAbsTestStatistic[_nGof];
529  _gofSplitMode.resize(_nGof);
530 
531  // Create array of regular fit contexts, containing subset of data and single fitCat PDF
532  catIter->Reset();
533  while ((type = (RooCatType*) catIter->Next())) {
534  // Retrieve the PDF for this simCat state
535  RooAbsPdf* pdf = simpdf->getPdf(type->GetName());
536  RooAbsData* dset = (RooAbsData*) dsetList->FindObject(type->GetName());
537 
538  if (pdf && dset && (0. != dset->sumEntries() || processEmptyDataSets())) {
539  ccoutI(Fitting) << "RooAbsTestStatistic::initSimMode: creating slave calculator #" << n << " for state " << type->GetName()
540  << " (" << dset->numEntries() << " dataset entries)" << endl;
541 
542 
543  // *** START HERE
544  // WVE HACK determine if we have a RooRealSumPdf and then treat it like a binned likelihood
545  RooAbsPdf* binnedPdf = 0 ;
546  Bool_t binnedL = kFALSE ;
547  if (pdf->getAttribute("BinnedLikelihood") && pdf->IsA()->InheritsFrom(RooRealSumPdf::Class())) {
548  // Simplest case: top-level of component is a RRSP
549  binnedPdf = pdf ;
550  binnedL = kTRUE ;
551  } else if (pdf->IsA()->InheritsFrom(RooProdPdf::Class())) {
552  // Default case: top-level pdf is a product of RRSP and other pdfs
553  RooFIter iter = ((RooProdPdf*)pdf)->pdfList().fwdIterator() ;
554  RooAbsArg* component ;
555  while ((component = iter.next())) {
556  if (component->getAttribute("BinnedLikelihood") && component->IsA()->InheritsFrom(RooRealSumPdf::Class())) {
557  binnedPdf = (RooAbsPdf*) component ;
558  binnedL = kTRUE ;
559  }
560  if (component->getAttribute("MAIN_MEASUREMENT")) {
561  // not really a binned pdf, but this prevents a (potentially) long list of subsidiary measurements to be passed to the slave calculator
562  binnedPdf = (RooAbsPdf*) component ;
563  }
564  }
565  }
566  // WVE END HACK
567  // Below here directly pass binnedPdf instead of PROD(binnedPdf,constraints) as constraints are evaluated elsewhere anyway
568  // and omitting them reduces model complexity and associated handling/cloning times
569  if (_splitRange && rangeName) {
570  _gofArray[n] = create(type->GetName(),type->GetName(),(binnedPdf?*binnedPdf:*pdf),*dset,*projDeps,
571  Form("%s_%s",rangeName,type->GetName()),addCoefRangeName,_nCPU*(_mpinterl?-1:1),_mpinterl,_verbose,_splitRange,binnedL);
572  } else {
573  _gofArray[n] = create(type->GetName(),type->GetName(),(binnedPdf?*binnedPdf:*pdf),*dset,*projDeps,
574  rangeName,addCoefRangeName,_nCPU,_mpinterl,_verbose,_splitRange,binnedL);
575  }
576  _gofArray[n]->setSimCount(_nGof);
577  // *** END HERE
578 
579  // Fill per-component split mode with Bulk Partition for now so that Auto will map to bulk-splitting of all components
580  if (_mpinterl==RooFit::Hybrid) {
581  if (dset->numEntries()<10) {
582  //cout << "RAT::initSim("<< GetName() << ") MP mode is auto, setting split mode for component "<< n << " to SimComponents"<< endl ;
583  _gofSplitMode[n] = RooFit::SimComponents;
584  _gofArray[n]->_mpinterl = RooFit::SimComponents;
585  } else {
586  //cout << "RAT::initSim("<< GetName() << ") MP mode is auto, setting split mode for component "<< n << " to BulkPartition"<< endl ;
587  _gofSplitMode[n] = RooFit::BulkPartition;
588  _gofArray[n]->_mpinterl = RooFit::BulkPartition;
589  }
590  }
591 
592  // Servers may have been redirected between instantiation and (deferred) initialization
593 
594  RooArgSet *actualParams = binnedPdf ? binnedPdf->getParameters(dset) : pdf->getParameters(dset);
595  RooArgSet* selTargetParams = (RooArgSet*) _paramSet.selectCommon(*actualParams);
596 
597  _gofArray[n]->recursiveRedirectServers(*selTargetParams);
598 
599  delete selTargetParams;
600  delete actualParams;
601 
602  ++n;
603 
604  } else {
605  if ((!dset || (0. != dset->sumEntries() && !processEmptyDataSets())) && pdf) {
606  if (_verbose) {
607  ccoutD(Fitting) << "RooAbsTestStatistic::initSimMode: state " << type->GetName()
608  << " has no data entries, no slave calculator created" << endl;
609  }
610  }
611  }
612  }
613  coutI(Fitting) << "RooAbsTestStatistic::initSimMode: created " << n << " slave calculators." << endl;
614 
615  dsetList->Delete(); // delete the content.
616  delete dsetList;
617  delete catIter;
618 }
619 
620 
621 ////////////////////////////////////////////////////////////////////////////////
622 /// Change dataset that is used to given one. If cloneData is kTRUE, a clone of
623 /// in the input dataset is made. If the test statistic was constructed with
624 /// a range specification on the data, the cloneData argument is ignore and
625 /// the data is always cloned.
626 Bool_t RooAbsTestStatistic::setData(RooAbsData& indata, Bool_t cloneData)
627 {
628  // Trigger refresh of likelihood offsets
629  if (isOffsetting()) {
630  enableOffsetting(kFALSE);
631  enableOffsetting(kTRUE);
632  }
633 
634  switch(operMode()) {
635  case Slave:
636  // Delegate to implementation
637  return setDataSlave(indata, cloneData);
638  case SimMaster:
639  // Forward to slaves
640  // cout << "RATS::setData(" << GetName() << ") SimMaster, calling setDataSlave() on slave nodes" << endl;
641  if (indata.canSplitFast()) {
642  for (Int_t i = 0; i < _nGof; ++i) {
643  RooAbsData* compData = indata.getSimData(_gofArray[i]->GetName());
644  _gofArray[i]->setDataSlave(*compData, cloneData);
645  }
646  } else if (0 == indata.numEntries()) {
647  // For an unsplit empty dataset, simply assign empty dataset to each component
648  for (Int_t i = 0; i < _nGof; ++i) {
649  _gofArray[i]->setDataSlave(indata, cloneData);
650  }
651  } else {
652  const RooAbsCategoryLValue& indexCat = static_cast<RooSimultaneous*>(_func)->indexCat();
653  TList* dlist = indata.split(indexCat, kTRUE);
654  if (!dlist) {
655  coutF(DataHandling) << "Tried to split '" << indata.GetName() << "' into categories of '" << indexCat.GetName()
656  << "', but splitting failed. Input data:" << std::endl;
657  indata.Print("V");
658  throw std::runtime_error("Error when setting up test statistic: dataset couldn't be split into categories.");
659  }
660 
661  for (Int_t i = 0; i < _nGof; ++i) {
662  RooAbsData* compData = (RooAbsData*) dlist->FindObject(_gofArray[i]->GetName());
663  // cout << "component data for index " << _gofArray[i]->GetName() << " is " << compData << endl;
664  if (compData) {
665  _gofArray[i]->setDataSlave(*compData,kFALSE,kTRUE);
666  } else {
667  coutE(DataHandling) << "RooAbsTestStatistic::setData(" << GetName() << ") ERROR: Cannot find component data for state " << _gofArray[i]->GetName() << endl;
668  }
669  }
670  }
671  break;
672  case MPMaster:
673  // Not supported
674  coutF(DataHandling) << "RooAbsTestStatistic::setData(" << GetName() << ") FATAL: setData() is not supported in multi-processor mode" << endl;
675  throw std::runtime_error("RooAbsTestStatistic::setData is not supported in MPMaster mode");
676  break;
677  }
678 
679  return kTRUE;
680 }
681 
682 
683 
684 void RooAbsTestStatistic::enableOffsetting(Bool_t flag)
685 {
686  // Apply internal value offsetting to control numeric precision
687  if (!_init) {
688  const_cast<RooAbsTestStatistic*>(this)->initialize() ;
689  }
690 
691  switch(operMode()) {
692  case Slave:
693  _doOffset = flag ;
694  // Clear offset if feature is disabled to that it is recalculated next time it is enabled
695  if (!_doOffset) {
696  _offset = 0 ;
697  _offsetCarry = 0;
698  }
699  setValueDirty() ;
700  break ;
701  case SimMaster:
702  _doOffset = flag;
703  for (Int_t i = 0; i < _nGof; ++i) {
704  _gofArray[i]->enableOffsetting(flag);
705  }
706  break ;
707  case MPMaster:
708  _doOffset = flag;
709  for (Int_t i = 0; i < _nCPU; ++i) {
710  _mpfeArray[i]->enableOffsetting(flag);
711  }
712  break;
713  }
714 }
715 
716 
717 Double_t RooAbsTestStatistic::getCarry() const
718 { return _evalCarry; }