laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauSimFitMaster.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2013 - 2013.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <cstdlib>
16 #include <iostream>
17 
18 #include "TMatrixD.h"
19 #include "TMessage.h"
20 #include "TMonitor.h"
21 #include "TObjArray.h"
22 #include "TObjString.h"
23 #include "TServerSocket.h"
24 #include "TSocket.h"
25 #include "TSystem.h"
26 
27 #include "LauAbsFitter.hh"
28 #include "LauFitNtuple.hh"
29 #include "LauFitter.hh"
30 #include "LauFormulaPar.hh"
31 #include "LauParameter.hh"
32 #include "LauParamFixed.hh"
33 #include "LauSimFitMaster.hh"
34 
35 
37 
38 
39 LauSimFitMaster::LauSimFitMaster( UInt_t numSlaves, UInt_t port ) :
40  nSlaves_(numSlaves),
41  reqPort_(port),
42  nParams_(0),
43  nFreeParams_(0),
44  withinAsymErrorCalc_(kFALSE),
45  numberOKFits_(0),
46  numberBadFits_(0),
47  fitStatus_(0),
48  iExpt_(0),
49  NLL_(0.0),
50  socketMonitor_(0),
51  messageFromSlave_(0),
52  fitNtuple_(0)
53 {
54  messagesToSlaves_.resize( nSlaves_ );
55  for ( UInt_t iSlave(0); iSlave < nSlaves_; ++iSlave ) {
56  messagesToSlaves_[iSlave] = new TMessage();
57  }
58 }
59 
61 {
62  delete socketMonitor_; socketMonitor_ = 0;
63 
64  // Tell all slaves that they are finished and delete corresponding socket
65  TString msgStr("Finish");
66  TMessage message( kMESS_STRING );
67  message.WriteTString(msgStr);
68  for ( std::vector<TSocket*>::iterator iter = sSlaves_.begin(); iter != sSlaves_.end(); ++iter ) {
69  (*iter)->Send(message);
70  (*iter)->Close();
71  delete (*iter);
72  }
73  sSlaves_.clear();
74 
75  // Remove the components created to apply constraints to fit parameters
76  for (std::vector<LauAbsRValue*>::iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter){
77  if ( !(*iter)->isLValue() ){
78  delete (*iter);
79  (*iter) = 0;
80  }
81  }
82  conVars_.clear();
83 
84  // Remove all fit parameters
85  for ( std::vector<LauParameter*>::iterator iter = params_.begin(); iter != params_.end(); ++iter ) {
86  delete *iter;
87  }
88  params_.clear();
89 
90  for ( std::vector<Double_t*>::iterator iter = vectorPar_.begin(); iter != vectorPar_.end(); ++iter ) {
91  delete[] (*iter);
92  }
93  vectorPar_.clear();
94 
96 
97  for ( std::vector<TMessage*>::iterator iter = messagesToSlaves_.begin(); iter != messagesToSlaves_.end(); ++iter ) {
98  delete (*iter);
99  }
100  messagesToSlaves_.clear();
101 
102  delete fitNtuple_;
103 }
104 
106 {
107  if ( socketMonitor_ != 0 ) {
108  std::cerr << "ERROR in LauSimFitMaster::initSockets : Sockets already initialised." << std::endl;
109  return;
110  }
111 
112  //initialise socket connection, then accept a connection and return a full-duplex communication socket.
113  socketMonitor_ = new TMonitor();
114 
115  TServerSocket *ss = new TServerSocket( reqPort_, kFALSE );
116  UInt_t actual_port = ss->GetLocalPort();
117 
118  std::cout << "INFO in LauSimFitMaster::initSockets : Waiting for connection with " << nSlaves_ << " workers on port " << actual_port << std::endl;
119 
120  sSlaves_.resize(nSlaves_);
121  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
122  sSlaves_[iSlave] = ss->Accept();
123  std::cout << " : Added slave " << iSlave << std::endl;
124  }
125 
126  // tell the clients to start
127  std::cout << "INFO in LauSimFitMaster::initSockets : Initialising slaves" << std::endl;
128  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
129 
130  TMessage message( kMESS_ANY );
131  message.WriteUInt(iSlave);
132  message.WriteUInt(nSlaves_);
133 
134  sSlaves_[iSlave]->Send(message);
135 
136  socketMonitor_->Add(sSlaves_[iSlave]);
137  }
138  std::cout << " : Now start fit\n" << std::endl;
139 
140  ss->Close();
141  delete ss;
142 }
143 
144 /*
145  * OLD VERSION THAT JUST GETS THE NAMES - COULD HAVE A SERIES OF EXCHANGES TO GET THE NAMES, INIT VALUES, RANGES, ETC. INSTEAD OF PASSING PARAMETERS
146  * THIS INCREASES THE GENERALITY OF THE CODE, I.E. THERE IS NO NEED FOR THE SLAVES TO KNOW ANY LAURA++ CLASS BUT THIS ONE, BUT MAKES IT RATHER MORE DENSE
147  * FOR THE MOMENT I WILL STICK WITH THE METHOD OF PASSING LAUPARAMETER OBJECTS AROUND AND CONSIDER GOING BACK TO THIS GENERAL METHOD ONCE EVERYTHING IS WORKING
148  *
149 void LauSimFitMaster::getParametersFromSlavesFirstTime()
150 {
151  slaveIndices_.resize( nSlaves_ );
152 
153  TSocket* sActive(0);
154  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
155  // Send a message to the slave, requesting the list of parameter names
156  TString msgStr = "Parameter Names";
157  TMessage message( kMESS_STRING );
158  message.WriteTString( msgStr );
159  sSlaves_[iSlave]->Send(message);
160 
161  // Wait to receive the response and check that it has come from the slave we just requested from
162  sActive = socketMonitor_->Select();
163  if ( sActive != sSlaves_[iSlave] ) {
164  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Received message from a different slave than expected!" << std::endl;
165  gSystem->Exit(1);
166  }
167 
168  // Read the object and extract the parameter names
169  sSlaves_[iSlave]->Recv( messageFromSlave_ );
170  TObjArray * objarray = dynamic_cast<TObjArray*>( messageFromSlave_->ReadObject( messageFromSlave_->GetClass() ) );
171  if ( ! objarray ) {
172  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Error reading parameter names from slave" << std::endl;
173  gSystem->Exit(1);
174  }
175 
176  Int_t nPars = objarray->GetEntries();
177  for ( Int_t iPar(0); iPar < nPars; ++iPar ) {
178  TObjString* objstring = dynamic_cast<TObjString*>( (*objarray)[iPar] );
179  if ( ! objstring ) {
180  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Error reading parameter names from slave" << std::endl;
181  gSystem->Exit(1);
182  }
183  TString parname = objstring->GetString();
184 
185  std::map< TString, UInt_t >::iterator iter = parIndices_.find( parname );
186  if ( iter != parIndices_.end() ) {
187  UInt_t index = iter->second;
188  slaveIndices_[iSlave].push_back( index );
189  } else {
190  UInt_t index = parIndices_.size();
191  parIndices_.insert( std::make_pair( parname, index ) );
192  parNames_.insert( std::make_pair( index, parname ) );
193  slaveIndices_[iSlave].push_back( index );
194  }
195  }
196 
197  delete objarray; objarray = 0;
198  delete messageFromSlave_; messageFromSlave_ = 0;
199  }
200 
201  UInt_t nPars = parNames_.size();
202  parValues_.resize( nPars );
203 }
204 */
205 
207 {
208  if ( socketMonitor_ == 0 ) {
209  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlaves : Sockets not initialised." << std::endl;
210  return;
211  }
212 
213  if ( params_.empty() ) {
215 
216  // Add variables to Gaussian constrain to a list
217  this->addConParameters();
218  } else {
220  }
221 }
222 
224 {
225  TSocket* sActive(0);
226 
227  // Construct a message, requesting the list of parameter names
228  TString msgStr = "Send Parameters";
229  TMessage message( kMESS_STRING );
230  message.WriteTString( msgStr );
231 
232  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
233  // Send the message to the slave
234  sSlaves_[iSlave]->Send(message);
235 
236  // Wait to receive the response and check that it has come from the slave we just requested from
237  sActive = socketMonitor_->Select();
238  if ( sActive != sSlaves_[iSlave] ) {
239  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Received message from a different slave than expected!" << std::endl;
240  gSystem->Exit(1);
241  }
242 
243  // Read the object and extract the parameter names
244  sSlaves_[iSlave]->Recv( messageFromSlave_ );
245  TObjArray * objarray = dynamic_cast<TObjArray*>( messageFromSlave_->ReadObject( messageFromSlave_->GetClass() ) );
246  if ( ! objarray ) {
247  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Error reading parameter names from slave" << std::endl;
248  gSystem->Exit(1);
249  }
250 
251  // We want to auto-delete the supplied parameters since we only copy their values in this case
252  objarray->SetOwner(kTRUE);
253 
254  const UInt_t nPars = objarray->GetEntries();
255  if ( nPars != slaveIndices_[iSlave].size() ) {
256  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Unexpected number of parameters received from slave" << std::endl;
257  gSystem->Exit(1);
258  }
259 
260  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
261  LauParameter* parameter = dynamic_cast<LauParameter*>( (*objarray)[iPar] );
262  if ( ! parameter ) {
263  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Error reading parameter from slave" << std::endl;
264  gSystem->Exit(1);
265  }
266 
267  TString parname = parameter->name();
268  Double_t parvalue = parameter->initValue();
269 
270  std::map< TString, UInt_t >::iterator iter = parIndices_.find( parname );
271  if ( iter == parIndices_.end() ) {
272  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Unexpected parameter name received from slave" << std::endl;
273  gSystem->Exit(1);
274  }
275 
276  const UInt_t index = iter->second;
277  if ( slaveIndices_[iSlave][iPar] != index ) {
278  std::cerr << "ERROR in LauSimFitMaster::updateParametersFromSlaves : Unexpected parameter received from slave" << std::endl;
279  gSystem->Exit(1);
280  }
281 
282  params_[index]->initValue( parvalue );
283  parValues_[index] = parvalue;
284  vectorPar_[iSlave][iPar] = parvalue;
285  this->checkParameter( parameter, index );
286  }
287 
288  delete objarray; objarray = 0;
290  }
291 }
292 
294 {
295  slaveIndices_.resize( nSlaves_ );
296  slaveFreeIndices_.resize( nSlaves_ );
297  vectorPar_.resize( nSlaves_ );
298  vectorRes_.resize( nSlaves_ );
299 
300  TSocket* sActive(0);
301 
302  // Construct a message, requesting the list of parameter names
303  TString msgStr = "Send Parameters";
304  TMessage message( kMESS_STRING );
305  message.WriteTString( msgStr );
306 
307  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
308  // Send the message to the slave
309  sSlaves_[iSlave]->Send(message);
310 
311  // Wait to receive the response and check that it has come from the slave we just requested from
312  sActive = socketMonitor_->Select();
313  if ( sActive != sSlaves_[iSlave] ) {
314  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Received message from a different slave than expected!" << std::endl;
315  gSystem->Exit(1);
316  }
317 
318  // Read the object and extract the parameter names
319  sSlaves_[iSlave]->Recv( messageFromSlave_ );
320  TObjArray * objarray = dynamic_cast<TObjArray*>( messageFromSlave_->ReadObject( messageFromSlave_->GetClass() ) );
321  if ( ! objarray ) {
322  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Error reading parameters from slave" << std::endl;
323  gSystem->Exit(1);
324  }
325 
326  const UInt_t nPars = objarray->GetEntries();
327 
328  vectorPar_[iSlave] = new Double_t[nPars];
329 
330  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
331  LauParameter* parameter = dynamic_cast<LauParameter*>( (*objarray)[iPar] );
332  if ( ! parameter ) {
333  std::cerr << "ERROR in LauSimFitMaster::getParametersFromSlavesFirstTime : Error reading parameter from slave" << std::endl;
334  gSystem->Exit(1);
335  }
336 
337  TString parname = parameter->name();
338  Double_t parvalue = parameter->initValue();
339  Bool_t parfixed = parameter->fixed();
340 
341  std::map< TString, UInt_t >::iterator iter = parIndices_.find( parname );
342  if ( iter != parIndices_.end() ) {
343  UInt_t index = iter->second;
344  slaveIndices_[iSlave].push_back( index );
345  if ( ! parfixed ) {
346  slaveFreeIndices_[iSlave].push_back( index );
347  }
348  this->checkParameter( parameter, index );
349  } else {
350  UInt_t index = parIndices_.size();
351  parIndices_.insert( std::make_pair( parname, index ) );
352  parNames_.insert( std::make_pair( index, parname ) );
353  slaveIndices_[iSlave].push_back( index );
354  if ( ! parfixed ) {
355  slaveFreeIndices_[iSlave].push_back( index );
356  }
357  params_.push_back( parameter );
358  parValues_.push_back( parvalue );
359  }
360  vectorPar_[iSlave][iPar] = parvalue;
361  }
362 
363  delete objarray; objarray = 0;
365  }
366 
367  nParams_ = params_.size();
368 }
369 
371 {
372  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
373  const std::vector<UInt_t>& indices = slaveIndices_[iSlave];
374 
375  std::cout << "INFO in LauSimFitMaster::printParInfo : Slave " << iSlave << " has the following parameters:\n";
376  for ( std::vector<UInt_t>::const_iterator iter = indices.begin(); iter != indices.end(); ++iter ) {
377  const TString& parName = parNames_.find(*iter)->second;
378  Double_t parValue = parValues_[*iter];
379  const LauParameter* par = params_[*iter];
380  if ( par->name() != parName || par->initValue() != parValue ) {
381  std::cerr << "ERROR in LauSimFitMaster::printParInfo : Discrepancy in parameter name and value records, this is very strange!!" << std::endl;
382  }
383 
384  std::cout << " : " << parName << " = " << parValue << " and has index " << *iter << "\n";
385  }
386 
387  std::cout << std::endl;
388  }
389 
390  std::cout << "INFO in LauSimFitMaster::printParInfo : " << "There are " << nParams_ << " parameters in total" << std::endl;
391 }
392 
393 void LauSimFitMaster::checkParameter( const LauParameter* param, UInt_t index ) const
394 {
395  const LauParameter* storedPar = params_[index];
396 
397  TString parName = storedPar->name();
398  if ( param->name() != parName ) {
399  std::cerr << "ERROR in LauSimFitMaster::checkParameter : Parameter name is different!! This shouldn't happen!!" << std::endl;
400  }
401  if ( param->initValue() != storedPar->initValue() ) {
402  std::cerr << "WARNING in LauSimFitMaster::checkParameter : Initial value for parameter " << parName << " is different, will use the value first set: " << storedPar->initValue() << std::endl;
403  }
404  if ( param->minValue() != storedPar->minValue() ) {
405  std::cerr << "WARNING in LauSimFitMaster::checkParameter : Minimum allowed value for parameter " << parName << " is different, will use the value first set: " << storedPar->minValue() << std::endl;
406  }
407  if ( param->maxValue() != storedPar->maxValue() ) {
408  std::cerr << "WARNING in LauSimFitMaster::checkParameter : Maximum allowed value for parameter " << parName << " is different, will use the value first set: " << storedPar->maxValue() << std::endl;
409  }
410  if ( param->fixed() != storedPar->fixed() ) {
411  std::cerr << "WARNING in LauSimFitMaster::checkParameter : Fixed/floating property of parameter " << parName << " is different, will use the value first set: " << (storedPar->fixed() ? "fixed" : "floating") << std::endl;
412  }
413  if ( param->secondStage() != storedPar->secondStage() ) {
414  std::cerr << "WARNING in LauSimFitMaster::checkParameter : Second stage property of parameter " << parName << " is different, will use the value first set: " << (storedPar->secondStage() ? "true" : "false") << std::endl;
415  }
416 }
417 
419 {
420  this->initSockets();
421 }
422 
423 void LauSimFitMaster::runSimFit( const TString& fitNtupleFileName, UInt_t nExpt, UInt_t firstExpt, Bool_t useAsymmErrors, Bool_t twoStageFit )
424 {
425  // Routine to perform the total fit.
426 
427  // First, initialise
428  this->initialise();
429 
430  std::cout << "INFO in LauSimFitMaster::runSimFit : First experiment = " << firstExpt << std::endl;
431  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of experiments = " << nExpt << std::endl;
432 
433  // Start the cumulative timer
434  cumulTimer_.Start();
435 
436  numberOKFits_ = 0, numberBadFits_ = 0;
437  fitStatus_ = -1;
438 
439  // Create and setup the fit results ntuple
440  std::cout << "INFO in LauSimFitMaster::runSimFit : Creating fit ntuple." << std::endl;
441  if (fitNtuple_ != 0) {delete fitNtuple_; fitNtuple_ = 0;}
442  fitNtuple_ = new LauFitNtuple(fitNtupleFileName, useAsymmErrors);
443 
444  // Loop over the number of experiments
445  for (iExpt_ = firstExpt; iExpt_ < (firstExpt+nExpt); ++iExpt_) {
446 
447  // Start the timer to see how long each fit takes
448  timer_.Start();
449 
450  // Instruct the slaves to read the data for this experiment
451  Bool_t readOK = this->readData();
452  if ( ! readOK ) {
453  std::cerr << "ERROR in LauSimFitMaster::runSimFit : One or more slaves reported problems with reading data for experiment " << iExpt_ << ", skipping..." << std::endl;
454  timer_.Stop();
455  continue;
456  }
457 
458  // Instruct the slaves to perform the caching
459  this->cacheInputData();
460 
461  // Do the fit
462  this->fitExpt( useAsymmErrors, twoStageFit );
463 
464  // Stop the timer and see how long the program took so far
465  timer_.Stop();
466  timer_.Print();
467 
468  // Instruct the slaves to finalise the results
469  this->finalise();
470 
471  // Keep track of how many fits succeeded or failed
472  if (fitStatus_ == 3) {
473  ++numberOKFits_;
474  } else {
475  ++numberBadFits_;
476  }
477  }
478 
479  // Print out total timing info.
480  std::cout << "INFO in LauSimFitMaster::runSimFit : Cumulative timing:" << std::endl;
481  cumulTimer_.Stop();
482  cumulTimer_.Print();
483 
484  // Print out stats on OK fits.
485  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of OK Fits = " << numberOKFits_ << std::endl;
486  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of Failed Fits = " << numberBadFits_ << std::endl;
487  Double_t fitEff(0.0);
488  if (nExpt != 0) {fitEff = numberOKFits_/(1.0*nExpt);}
489  std::cout << "INFO in LauSimFitMaster::runSimFit : Fit efficiency = " << fitEff*100.0 << "%." << std::endl;
490 
491  // Instruct the slaves to write out any fit results (ntuples etc...).
492  this->writeOutResults();
493 }
494 
496 {
497  if ( socketMonitor_ == 0 ) {
498  std::cerr << "ERROR in LauSimFitMaster::readData : Sockets not initialised." << std::endl;
499  return kFALSE;
500  }
501 
502  // Construct a message, requesting to read the data for the given experiment
503  TString msgStr("Read Expt");
504  TMessage message( kMESS_STRING );
505  message.WriteTString( msgStr );
506  message.WriteUInt( iExpt_ );
507 
508  // Send the message to the slaves
509  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
510  sSlaves_[iSlave]->Send(message);
511  }
512 
513  TSocket* sActive(0);
514  UInt_t responsesReceived(0);
515  Bool_t ok(kTRUE);
516  while ( responsesReceived != nSlaves_ ) {
517 
518  // Get the next queued response
519  sActive = socketMonitor_->Select();
520 
521  // Extract from the message the ID of the slave and the number of events read
522  sActive->Recv( messageFromSlave_ );
523  UInt_t iSlave(0);
524  UInt_t nEvents(0);
525  messageFromSlave_->ReadUInt( iSlave );
526  messageFromSlave_->ReadUInt( nEvents );
527 
528  if ( nEvents <= 0 ) {
529  std::cerr << "ERROR in LauSimFitMaster::readData : Slave " << iSlave << " reports no events found for experiment " << iExpt_ << std::endl;
530  ok = kFALSE;
531  } else {
532  std::cerr << "INFO in LauSimFitMaster::readData : Slave " << iSlave << " reports " << nEvents << " events found for experiment " << iExpt_ << std::endl;
533  }
534 
535  ++responsesReceived;
536  }
537 
538  return ok;
539 }
540 
542 {
543  if ( socketMonitor_ == 0 ) {
544  std::cerr << "ERROR in LauSimFitMaster::cacheInputData : Sockets not initialised." << std::endl;
545  return kFALSE;
546  }
547 
548  // Construct a message, requesting it to read the data for the given experiment
549  TString msgStr("Cache");
550  TMessage message( kMESS_STRING );
551  message.WriteTString( msgStr );
552 
553  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
554  // Send the message to the slave
555  sSlaves_[iSlave]->Send(message);
556  }
557 
558  TSocket* sActive(0);
559  UInt_t responsesReceived(0);
560  Bool_t allOK(kTRUE);
561  while ( responsesReceived != nSlaves_ ) {
562 
563  // Get the next queued response
564  sActive = socketMonitor_->Select();
565 
566  // Extract from the message the ID of the slave and the success/failure flag
567  sActive->Recv( messageFromSlave_ );
568  UInt_t iSlave(0);
569  Bool_t ok(kTRUE);
570  messageFromSlave_->ReadUInt( iSlave );
571  messageFromSlave_->ReadBool( ok );
572 
573  if ( ! ok ) {
574  std::cerr << "ERROR in LauSimFitMaster::cacheInputData : Slave " << iSlave << " reports an error performing caching" << std::endl;
575  allOK = kFALSE;
576  }
577 
578  ++responsesReceived;
579  }
580 
581  return allOK;
582 }
583 
585 {
586  this->getParametersFromSlaves();
587  this->printParInfo();
588 }
589 
590 void LauSimFitMaster::fitExpt( Bool_t useAsymmErrors, Bool_t twoStageFit )
591 {
592  // Routine to perform the actual fit for the given experiment
593 
594  // Instruct the salves to update initial fit parameters if required (e.g. if using random numbers).
595  this->checkInitFitParams();
596 
597  // Initialise the fitter
598  LauFitter::fitter()->useAsymmFitErrors( useAsymmErrors );
599  LauFitter::fitter()->twoStageFit( twoStageFit );
601 
604 
605  // Now ready for minimisation step
606  std::cout << "\nINFO in LauSimFitMaster::fitExpt : Start minimisation...\n";
607  std::pair<Int_t,Double_t> fitResult = LauFitter::fitter()->minimise();
608 
609  fitStatus_ = fitResult.first;
610  NLL_ = fitResult.second;
611 
612  // If we're doing a two stage fit we can now release (i.e. float)
613  // the 2nd stage parameters and re-fit
614  if ( twoStageFit ) {
615 
616  if ( fitStatus_ != 3 ) {
617  std::cerr << "ERROR in LauSimFitMaster:fitExpt : Not running second stage fit since first stage failed." << std::endl;
619  } else {
623  fitResult = LauFitter::fitter()->minimise();
624  }
625  }
626 
627  fitStatus_ = fitResult.first;
628  NLL_ = fitResult.second;
629  const TMatrixD& covMat = LauFitter::fitter()->covarianceMatrix();
630  covMatrix_.Clear();
631  covMatrix_.ResizeTo( covMat.GetNrows(), covMat.GetNcols() );
632  covMatrix_.SetMatrixArray( covMat.GetMatrixArray() );
633 
634  // Store the final fit results and errors into protected internal vectors that
635  // all sub-classes can use within their own finalFitResults implementation
636  // used below (e.g. putting them into an ntuple in a root file)
638 }
639 
640 void LauSimFitMaster::setParsFromMinuit(Double_t* par, Int_t npar)
641 {
642  // This function sets the internal parameters based on the values
643  // that Minuit is using when trying to minimise the total likelihood function.
644 
645  // MINOS reports different numbers of free parameters depending on the
646  // situation, so disable this check
647  if ( ! withinAsymErrorCalc_ ) {
648  if (static_cast<UInt_t>(npar) != nFreeParams_) {
649  std::cerr << "ERROR in LauSimFitMaster::setParsFromMinuit : Unexpected number of free parameters: " << npar << ".\n";
650  std::cerr << " Expected: " << nFreeParams_ << ".\n" << std::endl;
651  gSystem->Exit(EXIT_FAILURE);
652  }
653  }
654 
655  // Despite npar being the number of free parameters
656  // the par array actually contains all the parameters,
657  // free and floating...
658  // Update all the parameters with their new values.
659  // Change the value in the array to be sent out to the slaves and the
660  // parameters themselves (so that constraints are correctly calculated)
661  for (UInt_t i(0); i<nParams_; ++i) {
662  if (!params_[i]->fixed()) {
663  parValues_[i] = par[i];
664  params_[i]->value(par[i]);
665  }
666  }
667 }
668 
670 {
671  if ( socketMonitor_ == 0 ) {
672  std::cerr << "ERROR in LauSimFitMaster::getTotNegLogLikelihood : Sockets not initialised." << std::endl;
673  return 0.0;
674  }
675 
676  // Send current values of the parameters to the slaves.
677  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
678 
679  std::vector<UInt_t>& indices = slaveIndices_[iSlave];
680  std::vector<UInt_t>& freeIndices = slaveFreeIndices_[iSlave];
681  UInt_t nPars = indices.size();
682  UInt_t nFreePars = freeIndices.size();
683  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
684  vectorPar_[iSlave][iPar] = parValues_[ indices[iPar] ];
685  }
686 
687  TMessage* message = messagesToSlaves_[iSlave];
688  message->Reset( kMESS_ANY );
689  message->WriteUInt( nPars );
690  message->WriteUInt( nFreePars );
691  message->WriteFastArray( vectorPar_[iSlave], nPars );
692 
693  sSlaves_[iSlave]->Send(*message);
694  }
695 
696  Double_t negLogLike(0.0);
697  TSocket *sActive(0);
698  UInt_t responsesReceived(0);
699  while ( responsesReceived != nSlaves_ ) {
700 
701  sActive = socketMonitor_->Select();
702  sActive->Recv(messageFromSlave_);
703 
704  messageFromSlave_->ReadDouble( vectorRes_[responsesReceived] );
705 
706  negLogLike += vectorRes_[responsesReceived];
707 
708  ++responsesReceived;
709  }
710 
711  // Calculate any penalty terms from Gaussian constrained variables
712  if ( ! conVars_.empty() ){
713  negLogLike += this->getLogLikelihoodPenalty();
714  }
715 
716  return negLogLike;
717 }
718 
720 {
721  Double_t penalty(0.0);
722 
723  for ( std::vector<LauAbsRValue*>::const_iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter ) {
724  Double_t val = (*iter)->unblindValue();
725  Double_t mean = (*iter)->constraintMean();
726  Double_t width = (*iter)->constraintWidth();
727 
728  Double_t term = ( val - mean )*( val - mean );
729  penalty += term/( 2*width*width );
730  }
731 
732  return penalty;
733 }
734 
735 void LauSimFitMaster::addConstraint(const TString& formula, const std::vector<TString>& pars, const Double_t mean, const Double_t width)
736 {
737  StoreConstraints newCon;
738  newCon.formula_ = formula;
739  newCon.conPars_ = pars;
740  newCon.mean_ = mean;
741  newCon.width_ = width;
742  storeCon_.push_back(newCon);
743 }
744 
746 {
747  // Add penalties from the constraints to fit parameters
748 
749  // First, constraints on the fit parameters themselves
750  for ( std::vector<LauParameter*>::const_iterator iter = params_.begin(); iter != params_.end(); ++iter ) {
751  if ( (*iter)->gaussConstraint() ) {
752  conVars_.push_back( *iter );
753  std::cout << "INFO in LauSimFitMaster::addConParameters : Added Gaussian constraint to parameter "<< (*iter)->name() << std::endl;
754  }
755  }
756 
757  // Second, constraints on arbitrary combinations
758  for ( std::vector<StoreConstraints>::iterator iter = storeCon_.begin(); iter != storeCon_.end(); ++iter ) {
759  std::vector<TString> names = (*iter).conPars_;
760  std::vector<LauParameter*> params;
761  for ( std::vector<TString>::iterator iternames = names.begin(); iternames != names.end(); ++iternames ) {
762  for ( std::vector<LauParameter*>::const_iterator iterfit = params_.begin(); iterfit != params_.end(); ++iterfit ) {
763  if ( (*iternames) == (*iterfit)->name() ){
764  params.push_back(*iterfit);
765  }
766  }
767  }
768 
769  // If the parameters are not found, skip it
770  if ( params.size() != (*iter).conPars_.size() ) {
771  std::cerr << "WARNING in LauSimFitMaster::addConParameters: Could not find parameters to constrain in the formula... skipping" << std::endl;
772  continue;
773  }
774 
775  LauFormulaPar* formPar = new LauFormulaPar( (*iter).formula_, (*iter).formula_, params );
776  formPar->addGaussianConstraint( (*iter).mean_, (*iter).width_ );
777  conVars_.push_back(formPar);
778 
779  std::cout << "INFO in LauSimFitMaster::addConParameters : Added Gaussian constraint to formula\n";
780  std::cout << " : Formula: " << (*iter).formula_ << std::endl;
781  for ( std::vector<LauParameter*>::iterator iterparam = params.begin(); iterparam != params.end(); ++iterparam ) {
782  std::cout << " : Parameter: " << (*iterparam)->name() << std::endl;
783  }
784  }
785 
786 }
787 
789 {
790  if ( socketMonitor_ == 0 ) {
791  std::cerr << "ERROR in LauSimFitMaster::finalise : Sockets not initialised." << std::endl;
792  return kFALSE;
793  }
794 
795  // Prepare the covariance matrices
796  covMatrices_.resize( nSlaves_ );
797 
798  LauParamFixed pred;
799 
800  std::map<UInt_t,UInt_t> freeParIndices;
801 
802  UInt_t counter(0);
803  for ( UInt_t iPar(0); iPar < nParams_; ++iPar ) {
804  const LauParameter* par = params_[iPar];
805  if ( ! pred(par) ) {
806  freeParIndices.insert( std::make_pair(iPar,counter) );
807  ++counter;
808  }
809  }
810 
811  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
812  const UInt_t nPar = slaveIndices_[iSlave].size();
813 
814  std::vector<UInt_t> freeIndices;
815  freeIndices.reserve( nPar );
816 
817  for ( UInt_t iPar(0); iPar < nPar; ++iPar ) {
818  UInt_t index = slaveIndices_[iSlave][iPar];
819  std::map<UInt_t,UInt_t>::iterator freeIter = freeParIndices.find(index);
820  if ( freeIter == freeParIndices.end() ) {
821  continue;
822  }
823  UInt_t freeIndex = freeIter->second;
824  freeIndices.push_back( freeIndex );
825  }
826 
827  const UInt_t nFreePars = freeIndices.size();
828  TMatrixD& covMat = covMatrices_[iSlave];
829  covMat.ResizeTo( nFreePars, nFreePars );
830 
831  for ( UInt_t iPar(0); iPar < nFreePars; ++iPar ) {
832  for ( UInt_t jPar(0); jPar < nFreePars; ++jPar ) {
833  UInt_t i = freeIndices[iPar];
834  UInt_t j = freeIndices[jPar];
835  covMat( iPar, jPar ) = covMatrix_( i, j );
836  }
837  }
838  }
839 
840  // The array to hold the parameters
841  TObjArray array;
842 
843  // Send messages to all slaves containing the final parameters and fit status, NLL
844  // TODO - at present we lose the information on the correlations between the parameters that are unique to each slave
845  // - so should we store the full correlation matrix in an ntuple? along with all the parameters?
846  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
847 
848  array.Clear();
849 
850  std::vector<UInt_t>& indices = slaveIndices_[iSlave];
851  UInt_t nPars = indices.size();
852  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
853  array.Add( params_[ indices[iPar] ] );
854  }
855 
856  TMatrixD& covMat = covMatrices_[iSlave];
857 
858  TMessage* message = messagesToSlaves_[iSlave];
859  message->Reset( kMESS_OBJECT );
860  message->WriteInt( fitStatus_ );
861  message->WriteDouble( NLL_ );
862  message->WriteObject( &array );
863  message->WriteObject( &covMat );
864 
865  sSlaves_[iSlave]->Send(*message);
866  }
867 
868  TSocket *sActive(0);
869  UInt_t responsesReceived(0);
870  Bool_t allOK(kTRUE);
871  while ( responsesReceived != nSlaves_ ) {
872 
873  // Get the next queued response
874  sActive = socketMonitor_->Select();
875 
876  // Extract from the message the ID of the slave and the number of events read
877  sActive->Recv( messageFromSlave_ );
878  UInt_t iSlave(0);
879  Bool_t ok(kTRUE);
880  messageFromSlave_->ReadUInt( iSlave );
881  messageFromSlave_->ReadBool( ok );
882 
883  if ( ok ) {
884  TObjArray * objarray = dynamic_cast<TObjArray*>( messageFromSlave_->ReadObject( messageFromSlave_->GetClass() ) );
885  if ( ! objarray ) {
886  std::cerr << "ERROR in LauSimFitMaster::finalise : Error reading finalised parameters from slave" << std::endl;
887  allOK = kFALSE;
888  } else {
889  // We want to auto-delete the supplied parameters since we only copy their values in this case
890  objarray->SetOwner(kTRUE);
891 
892  const UInt_t nPars = objarray->GetEntries();
893  if ( nPars != slaveIndices_[iSlave].size() ) {
894  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected number of finalised parameters received from slave" << std::endl;
895  allOK = kFALSE;
896  } else {
897  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
898  LauParameter* parameter = dynamic_cast<LauParameter*>( (*objarray)[iPar] );
899  if ( ! parameter ) {
900  std::cerr << "ERROR in LauSimFitMaster::finalise : Error reading parameter from slave" << std::endl;
901  allOK = kFALSE;
902  continue;
903  }
904 
905  TString parname = parameter->name();
906 
907  std::map< TString, UInt_t >::iterator iter = parIndices_.find( parname );
908  if ( iter == parIndices_.end() ) {
909  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected parameter name received from slave" << std::endl;
910  allOK = kFALSE;
911  continue;
912  }
913 
914  const UInt_t index = iter->second;
915  if ( slaveIndices_[iSlave][iPar] != index ) {
916  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected parameter received from slave" << std::endl;
917  allOK = kFALSE;
918  continue;
919  }
920 
921  Double_t parvalue = parameter->value();
922  params_[index]->value( parvalue );
923  parValues_[index] = parvalue;
924  vectorPar_[iSlave][iPar] = parvalue;
925  }
926  }
927  delete objarray;
928  }
929  } else {
930  std::cerr << "ERROR in LauSimFitMaster::finalise : Slave " << iSlave << " reports an error performing finalisation" << std::endl;
931  allOK = kFALSE;
932  }
933 
934  ++responsesReceived;
935  }
936 
937  // Fill our ntuple as well
938  if ( fitNtuple_ != 0 ) {
939  for ( std::vector<LauParameter*>::iterator iter = params_.begin(); iter != params_.end(); ++iter ) {
940  if (!(*iter)->fixed()) {
941  (*iter)->updatePull();
942  }
943  }
944  std::vector<LauParameter> extraVars;
945  fitNtuple_->storeParsAndErrors( params_, extraVars );
948  }
949 
950  return allOK;
951 }
952 
954 {
955  if ( socketMonitor_ == 0 ) {
956  std::cerr << "ERROR in LauSimFitMaster::writeOutResults : Sockets not initialised." << std::endl;
957  return kFALSE;
958  }
959 
960  // Construct a message, requesting to write out the fit results
961  TString msgStr("Write Results");
962  TMessage message( kMESS_STRING );
963  message.WriteTString( msgStr );
964 
965  // Send the message to the slaves
966  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
967  sSlaves_[iSlave]->Send(message);
968  }
969 
970  TSocket *sActive(0);
971  UInt_t responsesReceived(0);
972  Bool_t allOK(kTRUE);
973  while ( responsesReceived != nSlaves_ ) {
974 
975  // Get the next queued response
976  sActive = socketMonitor_->Select();
977 
978  // Extract from the message the ID of the slave and the number of events read
979  sActive->Recv( messageFromSlave_ );
980  UInt_t iSlave(0);
981  Bool_t ok(kTRUE);
982  messageFromSlave_->ReadUInt( iSlave );
983  messageFromSlave_->ReadBool( ok );
984 
985  if ( ! ok ) {
986  std::cerr << "ERROR in LauSimFitMaster::writeOutResults : Slave " << iSlave << " reports an error performing finalisation" << std::endl;
987  allOK = kFALSE;
988  }
989 
990  ++responsesReceived;
991  }
992 
993  // Write out our ntuple as well
994  if (fitNtuple_ != 0) {
996  }
997 
998  return allOK;
999 }
1000 
TMonitor * socketMonitor_
Parallel setup monitor.
std::vector< std::vector< UInt_t > > slaveIndices_
Lists of indices for each slave.
std::vector< LauAbsRValue * > conVars_
Gaussian constraints.
Bool_t readData()
Instruct the slaves to read the input data for the given experiment.
File containing declaration of LauFormulaPar class.
void getParametersFromSlaves()
Determine/update the parameter initial values from all slaves.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
virtual void updateParameters()=0
Update the values and errors of the parameters based on the fit minimum.
TStopwatch timer_
The fit timer.
std::map< UInt_t, TString > parNames_
Reverse map of index in the values vector to parameter names.
void addConParameters()
Add parameters to the list of Gaussian constrained parameters.
File containing declaration of LauParamFixed class.
LauFitNtuple * fitNtuple_
The fit results ntuple.
Double_t maxValue() const
The maximum value allowed for the parameter.
Bool_t writeOutResults()
Instruct the slaves to write out the fit results.
ClassImp(LauAbsCoeffSet)
std::vector< LauParameter * > params_
Parameters.
const TString & name() const
The parameter name.
Class for defining combinations of fit parameter objects.
virtual Bool_t twoStageFit() const =0
Determine whether the two-stage fit is enabled.
File containing declaration of LauAbsFitter class.
Double_t minValue() const
The minimum value allowed for the parameter.
void getParametersFromSlavesFirstTime()
Determine the parameter names and initial values from all slaves.
UInt_t nParams_
The number of fit parameters.
std::vector< TSocket * > sSlaves_
Sockets for each of the slaves.
const UInt_t reqPort_
The requested port.
std::vector< StoreConstraints > storeCon_
Store the constraints for fit parameters until initialisation is complete.
void updateFitNtuple()
Update the fit ntuple.
std::vector< Double_t > vectorRes_
Likelihood values returned from the slaves.
std::map< TString, UInt_t > parIndices_
Map of parameter names to index in the values vector.
void addGaussianConstraint(Double_t newGaussMean, Double_t newGaussWidth)
Add a Gaussian constraint (or modify an existing one)
virtual void initialise(LauFitObject *fitObj, const std::vector< LauParameter * > &parameters)=0
Initialise the fitter, setting the information on the parameters.
virtual std::pair< Int_t, Double_t > minimise()=0
Perform the minimisation of the fit function.
Bool_t withinAsymErrorCalc_
Flag to indicate if the asymmetric error calculation (e.g. MINOS) is currently running.
File containing declaration of LauFitter class.
void initSockets()
Initialise socket connections for the slaves.
UInt_t numberBadFits_
The number of fit failures.
void writeOutFitResults()
Write out fit results.
TMatrixD covMatrix_
The fit covariance matrix.
Double_t getLogLikelihoodPenalty()
Calculate the penalty terms to the log likelihood from Gaussian constraints.
std::vector< std::vector< UInt_t > > slaveFreeIndices_
Lists of indices of free parameters for each slave.
void checkInitFitParams()
Instruct the slaves to update the initial fit parameter values, if required.
virtual void addConstraint(const TString &formula, const std::vector< TString > &pars, const Double_t mean, const Double_t width)
Store constraint information for fit parameters.
File containing declaration of LauParameter class.
Predicate to allow counting of the number of fixed parameters.
std::vector< Double_t * > vectorPar_
Parameter values to send to the slaves.
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit...
void updateParametersFromSlaves()
Update and verify the parameter initial values from all slaves.
std::vector< TMatrixD > covMatrices_
The covariance sub-matrices for each slave.
virtual UInt_t nParameters() const =0
Get the total number of fit parameters.
std::vector< TMessage * > messagesToSlaves_
Messages to slaves.
virtual void releaseSecondStageParameters()=0
Release parameters marked as &quot;second stage&quot;.
void storeCorrMatrix(UInt_t iExpt, Double_t NLL, Int_t fitStatus, const TMatrixD &covMatrix)
Store the correlation matrix and other fit information.
Definition: LauFitNtuple.cc:61
void checkParameter(const LauParameter *param, UInt_t index) const
Check for compatibility between two same-named parameters, which should therefore be identical...
Class for defining the fit parameter objects.
Definition: LauParameter.hh:35
void fitExpt(Bool_t useAsymmErrors, Bool_t twoStageFit)
Perform the fit for the current experiment.
Int_t fitStatus_
The status of the current fit.
The master process for simultaneous/combined fits.
static LauAbsFitter * fitter()
Method that provides access to the singleton fitter.
Definition: LauFitter.cc:34
Bool_t finalise()
Return the final parameters to the slaves and instruct them to perform their finalisation.
TStopwatch cumulTimer_
The total fit timer.
virtual UInt_t nFreeParameters() const =0
Get the number of floating fit parameters.
virtual void setParsFromMinuit(Double_t *par, Int_t npar)
This function sets the parameter values from Minuit.
std::vector< Double_t > parValues_
Parameter values.
const UInt_t nSlaves_
The number of slaves.
File containing declaration of LauSimFitMaster class.
void storeParsAndErrors(const std::vector< LauParameter * > &fitVars, const std::vector< LauParameter > &extraVars)
Store parameters and their errors.
void runSimFit(const TString &fitNtupleFileName, UInt_t nExpt, UInt_t firstExpt=0, Bool_t useAsymmErrors=kFALSE, Bool_t twoStageFit=kFALSE)
Run the fit.
Double_t initValue() const
The initial value of the parameter.
UInt_t nFreeParams_
The number of free fit parameters.
TMessage * messageFromSlave_
Message from slaves to the master.
virtual Bool_t useAsymmFitErrors() const =0
Determine whether calculation of asymmetric errors is enabled.
Bool_t cacheInputData()
Instruct the slaves to perform the caching.
Double_t NLL_
The negative log-likelihood.
Double_t value() const
The value of the parameter.
void initialise()
Initialise.
friend class LauFitNtuple
LauFitNtuple is a friend class.
virtual Double_t getTotNegLogLikelihood()
Calculate the new value of the negative log likelihood.
UInt_t iExpt_
The experiment number of the current fit.
virtual ~LauSimFitMaster()
Destructor.
UInt_t numberOKFits_
The number of successful fits.
File containing declaration of LauFitNtuple class.
virtual const TMatrixD & covarianceMatrix() const =0
Retrieve the fit covariance matrix.
void printParInfo() const
Print information on the parameters.