laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.

The base class for any task process for simultaneous/combined fits. More...

#include <LauSimFitTask.hh>

Inheritance diagram for LauSimFitTask:
LauFitObject LauAbsFitModel LauRooFitTask LauCPFitModel LauSimpleFitModel

Public Member Functions

 LauSimFitTask ()
 Constructor.
 
virtual ~LauSimFitTask ()
 Destructor.
 
UInt_t nTasks () const
 Obtain the number of tasks.
 
UInt_t taskId () const
 Obtain the ID number of this task.
 
virtual void runTask (const TString &dataFileName, const TString &dataTreeName, const TString &histFileName, const TString &tableFileName="", const TString &addressCoordinator="localhost", const UInt_t portCoordinator=9090)
 Start the task process for simultaneous fitting. More...
 
virtual void initialise ()=0
 Initialise the fit model. More...
 
virtual void setParsFromMinuit (Double_t *par, Int_t npar)=0
 This function sets the parameter values from Minuit. More...
 
virtual Double_t getTotNegLogLikelihood ()=0
 Calculates the total negative log-likelihood.
 
- Public Member Functions inherited from LauFitObject
virtual ~LauFitObject ()=default
 Destructor.
 
void useAsymmFitErrors (Bool_t useAsymmErrors)
 Turn on or off the computation of asymmetric errors (e.g. MINOS routine in Minuit) More...
 
Bool_t useAsymmFitErrors () const
 Report whether or not calculation of asymmetric errors is enabled.
 
void twoStageFit (Bool_t doTwoStageFit)
 Turn on or off the two stage fit. More...
 
Bool_t twoStageFit () const
 Report whether the two-stage fit is enabled.
 
virtual void withinAsymErrorCalc (const Bool_t inAsymErrCalc)
 Mark that the fit is calculating asymmetric errors. More...
 
Bool_t withinAsymErrorCalc () const
 Query whether the fit is calculating the asymmetric errors. More...
 
void setNExpts (UInt_t nExperiments, UInt_t firstExperiment, Bool_t toyExpts)
 Set the number of experiments, the first experiment, and whether this is toy. More...
 
UInt_t eventsPerExpt () const
 Obtain the total number of events in the current experiment.
 
UInt_t nExpt () const
 Obtain the number of experiments.
 
UInt_t firstExpt () const
 Obtain the number of the first experiment.
 
UInt_t iExpt () const
 Obtain the number of the current experiment.
 
Bool_t toyExpts () const
 Obtain whether this is toy.
 
void addConstraint (const TString &formula, const std::vector< TString > &pars, const Double_t mean, const Double_t width)
 Store constraint information for fit parameters. More...
 
void addFormulaConstraint (const TString &formula, const std::vector< TString > &pars, const Double_t mean, const Double_t width)
 Store constraint information for fit parameters. More...
 
void addMultiDimConstraint (const std::vector< TString > &pars, const TVectorD &means, const TMatrixD &covMat)
 Store n-dimensional constraint information for fit parameters. More...
 

Protected Member Functions

const LauFitNtuplefitNtuple () const
 Const access to the fit ntuple.
 
LauFitNtuplefitNtuple ()
 Access to the fit ntuple.
 
void connectToCoordinator (const TString &addressCoordinator, const UInt_t portCoordinator)
 Establish the connection to the coordinator process. More...
 
void processCoordinatorRequests ()
 Listen for requests from the coordinator and act accordingly.
 
virtual void setupResultsOutputs (const TString &histFileName, const TString &tableFileName)
 Setup saving of fit results to ntuple/LaTeX table etc. More...
 
virtual void prepareInitialParArray (TObjArray &array)=0
 Package the initial fit parameters for transmission to the coordinator. More...
 
virtual void finaliseExperiment (const LauAbsFitter::FitStatus &fitStat, const TObjArray *parsFromCoordinator, const TMatrixD *covMat, TObjArray &parsToCoordinator)=0
 Perform all finalisation actions. More...
 
virtual Bool_t verifyFitData (const TString &dataFileName, const TString &dataTreeName)=0
 Open the input file and verify that all required variables are present. More...
 
virtual UInt_t readExperimentData ()=0
 Read in the data for the current experiment. More...
 
virtual void cacheInputFitVars ()=0
 Cache the input data values to calculate the likelihood during the fit.
 
virtual void writeOutAllFitResults ()
 Write out any fit results.
 
- Protected Member Functions inherited from LauFitObject
 LauFitObject ()
 Constructor.
 
Bool_t checkRepetition (const std::vector< TString > &names, const ConstraintType conType)
 Check if parameters names for constraints have already been used elsewhere. More...
 
void generateConstraintMeans (std::vector< LauAbsRValue * > &conVars)
 Generate per-experiment mean for each Gaussian constraint. More...
 
const std::vector< FormulaConstraint > & formulaConstraints () const
 Const access to the formula constraints store.
 
std::vector< FormulaConstraint > & formulaConstraints ()
 Access to the formula constraints store.
 
const std::vector< MultiDimConstraint > & multiDimConstraints () const
 Const access to the ND constraints store.
 
std::vector< MultiDimConstraint > & multiDimConstraints ()
 Access to the ND constraints store.
 
const std::set< TString > & formulaConstrainedPars () const
 Const access to the parameter names used in formula constraints.
 
std::set< TString > & formulaConstrainedPars ()
 Access to the parameter names used in formula constraints.
 
const std::set< TString > & multiDimConstrainedPars () const
 Const access to the parameter names used in ND constraints.
 
std::set< TString > & multiDimConstrainedPars ()
 Access to the parameter names used in ND constraints.
 
void resetFitCounters ()
 Reset the good/bad fit counters.
 
void setCurrentExperiment (const UInt_t curExpt)
 Set the ID of the current experiment. More...
 
void startNewFit (const UInt_t nPars, const UInt_t nFreePars)
 Indicate the start of a new fit. More...
 
void eventsPerExpt (UInt_t nEvents)
 Set the number of events in the current experiment.
 
Double_t worstLogLike () const
 Access the worst log likelihood found so far.
 
void worstLogLike (const Double_t newWorstLogLike)
 Set a new value for the worst log likelihood. More...
 
void storeFitStatus (const LauAbsFitter::FitStatus &status, const TMatrixD &covMatrix)
 Store fit status information. More...
 
UInt_t nTotParams () const
 Access the total number of fit parameters.
 
UInt_t nFreeParams () const
 Access the total number of fit parameters.
 
const LauAbsFitter::FitStatusfitStatus () const
 Access the fit status information.
 
Double_t nll () const
 Access the current NLL value.
 
Double_t edm () const
 Access the current EDM value.
 
Int_t statusCode () const
 Access the fit status code.
 
const TMatrixD & covarianceMatrix () const
 Access the fit covariance matrix.
 
UInt_t numberOKFits () const
 Access the number of successful fits.
 
UInt_t numberBadFits () const
 Access the number of failed fits.
 

Private Member Functions

 LauSimFitTask (const LauSimFitTask &rhs)
 Copy constructor (not implemented)
 
LauSimFitTaskoperator= (const LauSimFitTask &rhs)
 Copy assignment operator (not implemented)
 

Private Attributes

TSocket * socketCoordinator_
 A socket to enable parallel setup.
 
TMessage * messageFromCoordinator_
 Message from coordinator to the tasks.
 
UInt_t taskId_
 Task id number.
 
UInt_t nTasks_
 The total number of tasks.
 
Double_t * parValues_
 Parameter values array (for reading from the coordinator)
 
LauFitNtuplefitNtuple_
 The fit ntuple.
 

Additional Inherited Members

- Protected Types inherited from LauFitObject
enum  ConstraintType { ConstraintType::Formula, ConstraintType::MultDim }
 Enumeration of the different types of constraint. More...
 

Detailed Description

The base class for any task process for simultaneous/combined fits.

Implementation of the JFit method described in arXiv:1409.5080 [physics.data-an].

This class acts as the base class from which tasks should inherit. This allows any fitting framework to plug in to the JFit method.

Definition at line 50 of file LauSimFitTask.hh.

Member Function Documentation

◆ connectToCoordinator()

void LauSimFitTask::connectToCoordinator ( const TString &  addressCoordinator,
const UInt_t  portCoordinator 
)
protected

Establish the connection to the coordinator process.

Parameters
[in]addressCoordinatorthe hostname of the machine running the coordinator process
[in]portCoordinatorthe port number on which the coordinator process is listening

Definition at line 114 of file LauSimFitTask.cc.

◆ finaliseExperiment()

virtual void LauSimFitTask::finaliseExperiment ( const LauAbsFitter::FitStatus fitStat,
const TObjArray *  parsFromCoordinator,
const TMatrixD *  covMat,
TObjArray &  parsToCoordinator 
)
protectedpure virtual

Perform all finalisation actions.

  • Receive the results of the fit from the coordinator
  • Perform any finalisation routines
  • Package the finalised fit parameters for transmission back to the coordinator
Parameters
[in]fitStatthe status of the fit, e.g. status code, EDM, NLL
[in]parsFromCoordinatorthe parameters at the fit minimum
[in]covMatthe fit covariance matrix
[out]parsToCoordinatorthe array to be filled with the finalised LauParameter objects

Implemented in LauAbsFitModel, and LauRooFitTask.

◆ initialise()

virtual void LauSimFitTask::initialise ( )
pure virtual

Initialise the fit model.

Each class that inherits from this one must implement this to do what is appropriate

Implemented in LauAbsFitModel, LauCPFitModel, LauSimpleFitModel, and LauRooFitTask.

◆ prepareInitialParArray()

virtual void LauSimFitTask::prepareInitialParArray ( TObjArray &  array)
protectedpure virtual

Package the initial fit parameters for transmission to the coordinator.

Parameters
[out]arraythe array to be filled with the LauParameter objects

Implemented in LauAbsFitModel, and LauRooFitTask.

◆ readExperimentData()

virtual UInt_t LauSimFitTask::readExperimentData ( )
protectedpure virtual

Read in the data for the current experiment.

Returns
the number of events read in

Implemented in LauAbsFitModel, and LauRooFitTask.

◆ runTask()

void LauSimFitTask::runTask ( const TString &  dataFileName,
const TString &  dataTreeName,
const TString &  histFileName,
const TString &  tableFileName = "",
const TString &  addressCoordinator = "localhost",
const UInt_t  portCoordinator = 9090 
)
virtual

Start the task process for simultaneous fitting.

Parameters
[in]dataFileNamethe name of the input data file
[in]dataTreeNamethe name of the tree containing the data
[in]histFileNamethe file name for the output histograms
[in]tableFileNamethe file name for the latex output file
[in]addressCoordinatorthe hostname of the machine running the coordinator process
[in]portCoordinatorthe port number on which the coordinator process is listening

Definition at line 61 of file LauSimFitTask.cc.

◆ setParsFromMinuit()

virtual void LauSimFitTask::setParsFromMinuit ( Double_t *  par,
Int_t  npar 
)
pure virtual

This function sets the parameter values from Minuit.

Parameters
[in]paran array storing the various parameter values
[in]nparthe number of free parameters

Implements LauFitObject.

Implemented in LauAbsFitModel, and LauRooFitTask.

◆ setupResultsOutputs()

void LauSimFitTask::setupResultsOutputs ( const TString &  histFileName,
const TString &  tableFileName 
)
protectedvirtual

Setup saving of fit results to ntuple/LaTeX table etc.

Provide here a default implementation that produces an ntuple only. Derived classes can override as they wish.

Parameters
[in]histFileNamethe file name for the output histograms
[in]tableFileNamethe file name for the latex output file

Reimplemented in LauAbsFitModel.

Definition at line 102 of file LauSimFitTask.cc.

◆ verifyFitData()

virtual Bool_t LauSimFitTask::verifyFitData ( const TString &  dataFileName,
const TString &  dataTreeName 
)
protectedpure virtual

Open the input file and verify that all required variables are present.

Parameters
[in]dataFileNamethe name of the input file
[in]dataTreeNamethe name of the input tree

Implemented in LauAbsFitModel, and LauRooFitTask.


The documentation for this class was generated from the following files: