laura is hosted by Hepforge, IPPP Durham
Laura++  v3r0p1
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->firstStage() != storedPar->firstStage() ) {
414  std::cerr << "WARNING in LauSimFitMaster::checkParameter : First stage property of parameter " << parName << " is different, will use the value first set: " << (storedPar->firstStage() ? "true" : "false") << std::endl;
415  }
416  if ( param->secondStage() != storedPar->secondStage() ) {
417  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;
418  }
419 }
420 
422 {
423  this->initSockets();
424 }
425 
426 void LauSimFitMaster::runSimFit( const TString& fitNtupleFileName, UInt_t nExpt, UInt_t firstExpt, Bool_t useAsymmErrors, Bool_t twoStageFit )
427 {
428  // Routine to perform the total fit.
429 
430  // First, initialise
431  this->initialise();
432 
433  std::cout << "INFO in LauSimFitMaster::runSimFit : First experiment = " << firstExpt << std::endl;
434  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of experiments = " << nExpt << std::endl;
435 
436  // Start the cumulative timer
437  cumulTimer_.Start();
438 
439  numberOKFits_ = 0, numberBadFits_ = 0;
440  fitStatus_ = -1;
441 
442  // Create and setup the fit results ntuple
443  std::cout << "INFO in LauSimFitMaster::runSimFit : Creating fit ntuple." << std::endl;
444  if (fitNtuple_ != 0) {delete fitNtuple_; fitNtuple_ = 0;}
445  fitNtuple_ = new LauFitNtuple(fitNtupleFileName, useAsymmErrors);
446 
447  // Loop over the number of experiments
448  for (iExpt_ = firstExpt; iExpt_ < (firstExpt+nExpt); ++iExpt_) {
449 
450  // Start the timer to see how long each fit takes
451  timer_.Start();
452 
453  // Instruct the slaves to read the data for this experiment
454  Bool_t readOK = this->readData();
455  if ( ! readOK ) {
456  std::cerr << "ERROR in LauSimFitMaster::runSimFit : One or more slaves reported problems with reading data for experiment " << iExpt_ << ", skipping..." << std::endl;
457  timer_.Stop();
458  continue;
459  }
460 
461  // Instruct the slaves to perform the caching
462  this->cacheInputData();
463 
464  // Do the fit
465  this->fitExpt( useAsymmErrors, twoStageFit );
466 
467  // Stop the timer and see how long the program took so far
468  timer_.Stop();
469  timer_.Print();
470 
471  // Instruct the slaves to finalise the results
472  this->finalise();
473 
474  // Keep track of how many fits succeeded or failed
475  if (fitStatus_ == 3) {
476  ++numberOKFits_;
477  } else {
478  ++numberBadFits_;
479  }
480  }
481 
482  // Print out total timing info.
483  std::cout << "INFO in LauSimFitMaster::runSimFit : Cumulative timing:" << std::endl;
484  cumulTimer_.Stop();
485  cumulTimer_.Print();
486 
487  // Print out stats on OK fits.
488  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of OK Fits = " << numberOKFits_ << std::endl;
489  std::cout << "INFO in LauSimFitMaster::runSimFit : Number of Failed Fits = " << numberBadFits_ << std::endl;
490  Double_t fitEff(0.0);
491  if (nExpt != 0) {fitEff = numberOKFits_/(1.0*nExpt);}
492  std::cout << "INFO in LauSimFitMaster::runSimFit : Fit efficiency = " << fitEff*100.0 << "%." << std::endl;
493 
494  // Instruct the slaves to write out any fit results (ntuples etc...).
495  this->writeOutResults();
496 }
497 
499 {
500  if ( socketMonitor_ == 0 ) {
501  std::cerr << "ERROR in LauSimFitMaster::readData : Sockets not initialised." << std::endl;
502  return kFALSE;
503  }
504 
505  // Construct a message, requesting to read the data for the given experiment
506  TString msgStr("Read Expt");
507  TMessage message( kMESS_STRING );
508  message.WriteTString( msgStr );
509  message.WriteUInt( iExpt_ );
510 
511  // Send the message to the slaves
512  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
513  sSlaves_[iSlave]->Send(message);
514  }
515 
516  TSocket* sActive(0);
517  UInt_t responsesReceived(0);
518  Bool_t ok(kTRUE);
519  while ( responsesReceived != nSlaves_ ) {
520 
521  // Get the next queued response
522  sActive = socketMonitor_->Select();
523 
524  // Extract from the message the ID of the slave and the number of events read
525  sActive->Recv( messageFromSlave_ );
526  UInt_t iSlave(0);
527  UInt_t nEvents(0);
528  messageFromSlave_->ReadUInt( iSlave );
529  messageFromSlave_->ReadUInt( nEvents );
530 
531  if ( nEvents <= 0 ) {
532  std::cerr << "ERROR in LauSimFitMaster::readData : Slave " << iSlave << " reports no events found for experiment " << iExpt_ << std::endl;
533  ok = kFALSE;
534  } else {
535  std::cerr << "INFO in LauSimFitMaster::readData : Slave " << iSlave << " reports " << nEvents << " events found for experiment " << iExpt_ << std::endl;
536  }
537 
538  ++responsesReceived;
539  }
540 
541  return ok;
542 }
543 
545 {
546  if ( socketMonitor_ == 0 ) {
547  std::cerr << "ERROR in LauSimFitMaster::cacheInputData : Sockets not initialised." << std::endl;
548  return kFALSE;
549  }
550 
551  // Construct a message, requesting it to read the data for the given experiment
552  TString msgStr("Cache");
553  TMessage message( kMESS_STRING );
554  message.WriteTString( msgStr );
555 
556  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
557  // Send the message to the slave
558  sSlaves_[iSlave]->Send(message);
559  }
560 
561  TSocket* sActive(0);
562  UInt_t responsesReceived(0);
563  Bool_t allOK(kTRUE);
564  while ( responsesReceived != nSlaves_ ) {
565 
566  // Get the next queued response
567  sActive = socketMonitor_->Select();
568 
569  // Extract from the message the ID of the slave and the success/failure flag
570  sActive->Recv( messageFromSlave_ );
571  UInt_t iSlave(0);
572  Bool_t ok(kTRUE);
573  messageFromSlave_->ReadUInt( iSlave );
574  messageFromSlave_->ReadBool( ok );
575 
576  if ( ! ok ) {
577  std::cerr << "ERROR in LauSimFitMaster::cacheInputData : Slave " << iSlave << " reports an error performing caching" << std::endl;
578  allOK = kFALSE;
579  }
580 
581  ++responsesReceived;
582  }
583 
584  return allOK;
585 }
586 
588 {
589  this->getParametersFromSlaves();
590  this->printParInfo();
591 }
592 
593 void LauSimFitMaster::fitExpt( Bool_t useAsymmErrors, Bool_t twoStageFit )
594 {
595  // Routine to perform the actual fit for the given experiment
596 
597  // Instruct the salves to update initial fit parameters if required (e.g. if using random numbers).
598  this->checkInitFitParams();
599 
600  // Initialise the fitter
601  LauFitter::fitter()->useAsymmFitErrors( useAsymmErrors );
602  LauFitter::fitter()->twoStageFit( twoStageFit );
604 
607 
608  // Now ready for minimisation step
609  std::cout << "\nINFO in LauSimFitMaster::fitExpt : Start minimisation...\n";
610  std::pair<Int_t,Double_t> fitResult = LauFitter::fitter()->minimise();
611 
612  fitStatus_ = fitResult.first;
613  NLL_ = fitResult.second;
614 
615  // If we're doing a two stage fit we can now release (i.e. float)
616  // the 2nd stage parameters and re-fit
617  if ( twoStageFit ) {
618 
619  if ( fitStatus_ != 3 ) {
620  std::cerr << "ERROR in LauSimFitMaster:fitExpt : Not running second stage fit since first stage failed." << std::endl;
622  } else {
627  fitResult = LauFitter::fitter()->minimise();
628  }
629  }
630 
631  fitStatus_ = fitResult.first;
632  NLL_ = fitResult.second;
633  const TMatrixD& covMat = LauFitter::fitter()->covarianceMatrix();
634  covMatrix_.Clear();
635  covMatrix_.ResizeTo( covMat.GetNrows(), covMat.GetNcols() );
636  covMatrix_.SetMatrixArray( covMat.GetMatrixArray() );
637 
638  // Store the final fit results and errors into protected internal vectors that
639  // all sub-classes can use within their own finalFitResults implementation
640  // used below (e.g. putting them into an ntuple in a root file)
643 }
644 
645 void LauSimFitMaster::setParsFromMinuit(Double_t* par, Int_t npar)
646 {
647  // This function sets the internal parameters based on the values
648  // that Minuit is using when trying to minimise the total likelihood function.
649 
650  // MINOS reports different numbers of free parameters depending on the
651  // situation, so disable this check
652  if ( ! withinAsymErrorCalc_ ) {
653  if (static_cast<UInt_t>(npar) != nFreeParams_) {
654  std::cerr << "ERROR in LauSimFitMaster::setParsFromMinuit : Unexpected number of free parameters: " << npar << ".\n";
655  std::cerr << " Expected: " << nFreeParams_ << ".\n" << std::endl;
656  gSystem->Exit(EXIT_FAILURE);
657  }
658  }
659 
660  // Despite npar being the number of free parameters
661  // the par array actually contains all the parameters,
662  // free and floating...
663  // Update all the parameters with their new values.
664  // Change the value in the array to be sent out to the slaves and the
665  // parameters themselves (so that constraints are correctly calculated)
666  for (UInt_t i(0); i<nParams_; ++i) {
667  if (!params_[i]->fixed()) {
668  parValues_[i] = par[i];
669  params_[i]->value(par[i]);
670  }
671  }
672 }
673 
675 {
676  if ( socketMonitor_ == 0 ) {
677  std::cerr << "ERROR in LauSimFitMaster::getTotNegLogLikelihood : Sockets not initialised." << std::endl;
678  return 0.0;
679  }
680 
681  // Send current values of the parameters to the slaves.
682  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
683 
684  std::vector<UInt_t>& indices = slaveIndices_[iSlave];
685  std::vector<UInt_t>& freeIndices = slaveFreeIndices_[iSlave];
686  UInt_t nPars = indices.size();
687  UInt_t nFreePars = freeIndices.size();
688  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
689  vectorPar_[iSlave][iPar] = parValues_[ indices[iPar] ];
690  }
691 
692  TMessage* message = messagesToSlaves_[iSlave];
693  message->Reset( kMESS_ANY );
694  message->WriteUInt( nPars );
695  message->WriteUInt( nFreePars );
696  message->WriteFastArray( vectorPar_[iSlave], nPars );
697 
698  sSlaves_[iSlave]->Send(*message);
699  }
700 
701  Double_t negLogLike(0.0);
702  TSocket *sActive(0);
703  UInt_t responsesReceived(0);
704  while ( responsesReceived != nSlaves_ ) {
705 
706  sActive = socketMonitor_->Select();
707  sActive->Recv(messageFromSlave_);
708 
709  messageFromSlave_->ReadDouble( vectorRes_[responsesReceived] );
710 
711  negLogLike += vectorRes_[responsesReceived];
712 
713  ++responsesReceived;
714  }
715 
716  // Calculate any penalty terms from Gaussian constrained variables
717  if ( ! conVars_.empty() ){
718  negLogLike += this->getLogLikelihoodPenalty();
719  }
720 
721  return negLogLike;
722 }
723 
725 {
726  Double_t penalty(0.0);
727 
728  for ( std::vector<LauAbsRValue*>::const_iterator iter = conVars_.begin(); iter != conVars_.end(); ++iter ) {
729  Double_t val = (*iter)->value();
730  Double_t mean = (*iter)->constraintMean();
731  Double_t width = (*iter)->constraintWidth();
732 
733  Double_t term = ( val - mean )*( val - mean );
734  penalty += term/( 2*width*width );
735  }
736 
737  return penalty;
738 }
739 
740 void LauSimFitMaster::addConstraint(const TString& formula, const std::vector<TString>& pars, const Double_t mean, const Double_t width)
741 {
742  StoreConstraints newCon;
743  newCon.formula_ = formula;
744  newCon.conPars_ = pars;
745  newCon.mean_ = mean;
746  newCon.width_ = width;
747  storeCon_.push_back(newCon);
748 }
749 
751 {
752  // Add penalties from the constraints to fit parameters
753 
754  // First, constraints on the fit parameters themselves
755  for ( std::vector<LauParameter*>::const_iterator iter = params_.begin(); iter != params_.end(); ++iter ) {
756  if ( (*iter)->gaussConstraint() ) {
757  conVars_.push_back( *iter );
758  std::cout << "INFO in LauSimFitMaster::addConParameters : Added Gaussian constraint to parameter "<< (*iter)->name() << std::endl;
759  }
760  }
761 
762  // Second, constraints on arbitrary combinations
763  for ( std::vector<StoreConstraints>::iterator iter = storeCon_.begin(); iter != storeCon_.end(); ++iter ) {
764  std::vector<TString> names = (*iter).conPars_;
765  std::vector<LauParameter*> params;
766  for ( std::vector<TString>::iterator iternames = names.begin(); iternames != names.end(); ++iternames ) {
767  for ( std::vector<LauParameter*>::const_iterator iterfit = params_.begin(); iterfit != params_.end(); ++iterfit ) {
768  if ( (*iternames) == (*iterfit)->name() ){
769  params.push_back(*iterfit);
770  }
771  }
772  }
773 
774  // If the parameters are not found, skip it
775  if ( params.size() != (*iter).conPars_.size() ) {
776  std::cerr << "WARNING in LauSimFitMaster::addConParameters: Could not find parameters to constrain in the formula... skipping" << std::endl;
777  continue;
778  }
779 
780  LauFormulaPar* formPar = new LauFormulaPar( (*iter).formula_, (*iter).formula_, params );
781  formPar->addGaussianConstraint( (*iter).mean_, (*iter).width_ );
782  conVars_.push_back(formPar);
783 
784  std::cout << "INFO in LauSimFitMaster::addConParameters : Added Gaussian constraint to formula\n";
785  std::cout << " : Formula: " << (*iter).formula_ << std::endl;
786  for ( std::vector<LauParameter*>::iterator iterparam = params.begin(); iterparam != params.end(); ++iterparam ) {
787  std::cout << " : Parameter: " << (*iterparam)->name() << std::endl;
788  }
789  }
790 
791 }
792 
794 {
795  if ( socketMonitor_ == 0 ) {
796  std::cerr << "ERROR in LauSimFitMaster::finalise : Sockets not initialised." << std::endl;
797  return kFALSE;
798  }
799 
800  // Prepare the covariance matrices
801  covMatrices_.resize( nSlaves_ );
802 
803  LauParamFixed pred;
804 
805  std::map<UInt_t,UInt_t> freeParIndices;
806 
807  UInt_t counter(0);
808  for ( UInt_t iPar(0); iPar < nParams_; ++iPar ) {
809  const LauParameter* par = params_[iPar];
810  if ( ! pred(par) ) {
811  freeParIndices.insert( std::make_pair(iPar,counter) );
812  ++counter;
813  }
814  }
815 
816  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
817  const UInt_t nPar = slaveIndices_[iSlave].size();
818 
819  std::vector<UInt_t> freeIndices;
820  freeIndices.reserve( nPar );
821 
822  for ( UInt_t iPar(0); iPar < nPar; ++iPar ) {
823  UInt_t index = slaveIndices_[iSlave][iPar];
824  std::map<UInt_t,UInt_t>::iterator freeIter = freeParIndices.find(index);
825  if ( freeIter == freeParIndices.end() ) {
826  continue;
827  }
828  UInt_t freeIndex = freeIter->second;
829  freeIndices.push_back( freeIndex );
830  }
831 
832  const UInt_t nFreePars = freeIndices.size();
833  TMatrixD& covMat = covMatrices_[iSlave];
834  covMat.ResizeTo( nFreePars, nFreePars );
835 
836  for ( UInt_t iPar(0); iPar < nFreePars; ++iPar ) {
837  for ( UInt_t jPar(0); jPar < nFreePars; ++jPar ) {
838  UInt_t i = freeIndices[iPar];
839  UInt_t j = freeIndices[jPar];
840  covMat( iPar, jPar ) = covMatrix_( i, j );
841  }
842  }
843  }
844 
845  // The array to hold the parameters
846  TObjArray array;
847 
848  // Send messages to all slaves containing the final parameters and fit status, NLL
849  // TODO - at present we lose the information on the correlations between the parameters that are unique to each slave
850  // - so should we store the full correlation matrix in an ntuple? along with all the parameters?
851  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
852 
853  array.Clear();
854 
855  std::vector<UInt_t>& indices = slaveIndices_[iSlave];
856  UInt_t nPars = indices.size();
857  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
858  array.Add( params_[ indices[iPar] ] );
859  }
860 
861  TMatrixD& covMat = covMatrices_[iSlave];
862 
863  TMessage* message = messagesToSlaves_[iSlave];
864  message->Reset( kMESS_OBJECT );
865  message->WriteInt( fitStatus_ );
866  message->WriteDouble( NLL_ );
867  message->WriteObject( &array );
868  message->WriteObject( &covMat );
869 
870  sSlaves_[iSlave]->Send(*message);
871  }
872 
873  TSocket *sActive(0);
874  UInt_t responsesReceived(0);
875  Bool_t allOK(kTRUE);
876  while ( responsesReceived != nSlaves_ ) {
877 
878  // Get the next queued response
879  sActive = socketMonitor_->Select();
880 
881  // Extract from the message the ID of the slave and the number of events read
882  sActive->Recv( messageFromSlave_ );
883  UInt_t iSlave(0);
884  Bool_t ok(kTRUE);
885  messageFromSlave_->ReadUInt( iSlave );
886  messageFromSlave_->ReadBool( ok );
887 
888  if ( ok ) {
889  TObjArray * objarray = dynamic_cast<TObjArray*>( messageFromSlave_->ReadObject( messageFromSlave_->GetClass() ) );
890  if ( ! objarray ) {
891  std::cerr << "ERROR in LauSimFitMaster::finalise : Error reading finalised parameters from slave" << std::endl;
892  allOK = kFALSE;
893  } else {
894  // We want to auto-delete the supplied parameters since we only copy their values in this case
895  objarray->SetOwner(kTRUE);
896 
897  const UInt_t nPars = objarray->GetEntries();
898  if ( nPars != slaveIndices_[iSlave].size() ) {
899  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected number of finalised parameters received from slave" << std::endl;
900  allOK = kFALSE;
901  } else {
902  for ( UInt_t iPar(0); iPar < nPars; ++iPar ) {
903  LauParameter* parameter = dynamic_cast<LauParameter*>( (*objarray)[iPar] );
904  if ( ! parameter ) {
905  std::cerr << "ERROR in LauSimFitMaster::finalise : Error reading parameter from slave" << std::endl;
906  allOK = kFALSE;
907  continue;
908  }
909 
910  TString parname = parameter->name();
911 
912  std::map< TString, UInt_t >::iterator iter = parIndices_.find( parname );
913  if ( iter == parIndices_.end() ) {
914  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected parameter name received from slave" << std::endl;
915  allOK = kFALSE;
916  continue;
917  }
918 
919  const UInt_t index = iter->second;
920  if ( slaveIndices_[iSlave][iPar] != index ) {
921  std::cerr << "ERROR in LauSimFitMaster::finalise : Unexpected parameter received from slave" << std::endl;
922  allOK = kFALSE;
923  continue;
924  }
925 
926  Double_t parvalue = parameter->value();
927  params_[index]->value( parvalue );
928  parValues_[index] = parvalue;
929  vectorPar_[iSlave][iPar] = parvalue;
930  }
931  }
932  delete objarray;
933  }
934  } else {
935  std::cerr << "ERROR in LauSimFitMaster::finalise : Slave " << iSlave << " reports an error performing finalisation" << std::endl;
936  allOK = kFALSE;
937  }
938 
939  ++responsesReceived;
940  }
941 
942  // Fill our ntuple as well
943  if ( fitNtuple_ != 0 ) {
944  for ( std::vector<LauParameter*>::iterator iter = params_.begin(); iter != params_.end(); ++iter ) {
945  if (!(*iter)->fixed()) {
946  (*iter)->updatePull();
947  }
948  }
949  std::vector<LauParameter> extraVars;
950  fitNtuple_->storeParsAndErrors( params_, extraVars );
953  }
954 
955  return allOK;
956 }
957 
959 {
960  if ( socketMonitor_ == 0 ) {
961  std::cerr << "ERROR in LauSimFitMaster::writeOutResults : Sockets not initialised." << std::endl;
962  return kFALSE;
963  }
964 
965  // Construct a message, requesting to write out the fit results
966  TString msgStr("Write Results");
967  TMessage message( kMESS_STRING );
968  message.WriteTString( msgStr );
969 
970  // Send the message to the slaves
971  for ( UInt_t iSlave(0); iSlave<nSlaves_; ++iSlave ) {
972  sSlaves_[iSlave]->Send(message);
973  }
974 
975  TSocket *sActive(0);
976  UInt_t responsesReceived(0);
977  Bool_t allOK(kTRUE);
978  while ( responsesReceived != nSlaves_ ) {
979 
980  // Get the next queued response
981  sActive = socketMonitor_->Select();
982 
983  // Extract from the message the ID of the slave and the number of events read
984  sActive->Recv( messageFromSlave_ );
985  UInt_t iSlave(0);
986  Bool_t ok(kTRUE);
987  messageFromSlave_->ReadUInt( iSlave );
988  messageFromSlave_->ReadBool( ok );
989 
990  if ( ! ok ) {
991  std::cerr << "ERROR in LauSimFitMaster::writeOutResults : Slave " << iSlave << " reports an error performing finalisation" << std::endl;
992  allOK = kFALSE;
993  }
994 
995  ++responsesReceived;
996  }
997 
998  // Write out our ntuple as well
999  if (fitNtuple_ != 0) {
1001  }
1002 
1003  return allOK;
1004 }
1005 
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.
Double_t mean_
The mean value of the Gaussian constraint to be applied.
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 LauFitter class.
void initSockets()
Initialise socket connections for the slaves.
Struct to store constraint information until the fit is run.
std::vector< TString > conPars_
The list of LauParameter names to be used in the LauFormulaPar.
UInt_t numberBadFits_
The number of fit failures.
Bool_t firstStage() const
Check whether the parameter should be floated only in the first stage of a two stage fit...
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.
virtual void fixFirstStageParameters()=0
Fix parameters marked as &quot;first stage&quot;.
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:64
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:34
void fitExpt(Bool_t useAsymmErrors, Bool_t twoStageFit)
Perform the fit for the current experiment.
TString formula_
The formula to be used in the LauFormulaPar.
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.
Double_t width_
The width of the Gaussian constraint to be applied.
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.
virtual void releaseFirstStageParameters()=0
Release parameters marked as &quot;first stage&quot;.
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.