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