laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauSimFitCoordinator.hh
Go to the documentation of this file.
1 
2 /*
3 Copyright 2013 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
37 #ifndef LAU_SIM_FIT_COORDINATOR
38 #define LAU_SIM_FIT_COORDINATOR
39 
40 #include "LauFitObject.hh"
41 
42 #include "TMatrixD.h"
43 #include "TStopwatch.h"
44 #include "TString.h"
45 
46 #include <map>
47 #include <vector>
48 
49 class TMessage;
50 class TMonitor;
51 class TSocket;
52 class LauAbsRValue;
53 class LauParameter;
54 class LauFitNtuple;
55 
57 
58  public:
60 
64  LauSimFitCoordinator( UInt_t numTasks, UInt_t port = 9090 );
65 
67  virtual ~LauSimFitCoordinator();
68 
70 
75  void runSimFit( const TString& fitNtupleFileName,
76  const Bool_t useAsymmErrors = kFALSE,
77  const Bool_t doTwoStageFit = kFALSE );
78 
80 
88  virtual void withinAsymErrorCalc( const Bool_t inAsymErrCalc );
89 
90  // Need to unshadow the query method defined in the base class
92 
94 
101  virtual void setParsFromMinuit( Double_t* par, Int_t npar );
102 
104 
108  virtual Double_t getTotNegLogLikelihood();
109 
110  protected:
112  void printParInfo() const;
113 
115  void initialise();
116 
118  void initSockets();
119 
121  void getParametersFromTasks();
122 
125 
128 
130  void checkParameter( const LauParameter* param, UInt_t index ) const;
131 
133  void addConParameters();
134 
136  Double_t getLogLikelihoodPenalty();
137 
139 
142  Bool_t readData();
143 
145 
148  Bool_t cacheInputData();
149 
151  void fitExpt();
152 
154  void checkInitFitParams();
155 
157  Bool_t finalise();
158 
160  Bool_t writeOutResults();
161 
162  private:
165 
168 
170  const UInt_t nTasks_;
171 
173  const UInt_t reqPort_;
174 
176  std::vector<TMatrixD> covMatrices_;
177 
179  TMonitor* socketMonitor_;
180 
182  std::vector<TSocket*> socketTasks_;
183 
185  std::vector<TMessage*> messagesToTasks_;
186 
188  TMessage* messageFromTask_;
189 
191  std::map<TString, UInt_t> parIndices_;
192 
194  std::map<UInt_t, TString> parNames_;
195 
197  std::vector<LauParameter*> params_;
198 
200  std::vector<LauAbsRValue*> conVars_;
201 
203  std::vector<Double_t> parValues_;
204 
206  std::vector<std::vector<UInt_t>> taskIndices_;
207 
209  std::vector<std::vector<UInt_t>> taskFreeIndices_;
210 
212  std::vector<Double_t*> vectorPar_;
213 
215  std::vector<Double_t> vectorRes_;
216 
218  TStopwatch timer_;
219 
221  TStopwatch cumulTimer_;
222 
225 
226  ClassDef( LauSimFitCoordinator, 0 );
227 };
228 
229 #endif
void fitExpt()
Perform the fit for the current experiment.
Bool_t finalise()
Return the final parameters to the tasks and instruct them to perform their finalisation.
virtual void setParsFromMinuit(Double_t *par, Int_t npar)
This function sets the parameter values from Minuit.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
std::vector< TMatrixD > covMatrices_
The covariance sub-matrices for each task.
std::vector< Double_t > parValues_
Parameter values.
void printParInfo() const
Print information on the parameters.
std::vector< std::vector< UInt_t > > taskFreeIndices_
Lists of indices of free parameters for each task.
const UInt_t reqPort_
The requested port.
TStopwatch timer_
The fit timer.
std::vector< LauParameter * > params_
Parameters.
Bool_t readData()
Instruct the tasks to read the input data for the given experiment.
LauSimFitCoordinator(const LauSimFitCoordinator &rhs)
Copy constructor (not implemented)
std::vector< std::vector< UInt_t > > taskIndices_
Lists of indices for each task.
void initSockets()
Initialise socket connections for the tasks.
void getParametersFromTasksFirstTime()
Determine the parameter names and initial values from all tasks.
std::vector< LauAbsRValue * > conVars_
Gaussian constraints.
Bool_t cacheInputData()
Instruct the tasks to perform the caching.
TMessage * messageFromTask_
Message from tasks to the coordinator.
Bool_t withinAsymErrorCalc() const
Query whether the fit is calculating the asymmetric errors.
Definition: LauFitObject.hh:94
void getParametersFromTasks()
Determine/update the parameter initial values from all tasks.
Class to store the results from the fit into an ntuple.
Definition: LauFitNtuple.hh:58
virtual ~LauSimFitCoordinator()
Destructor.
LauSimFitCoordinator(UInt_t numTasks, UInt_t port=9090)
Constructor.
void addConParameters()
Add parameters to the list of Gaussian constrained parameters.
std::vector< Double_t * > vectorPar_
Parameter values to send to the tasks.
void checkParameter(const LauParameter *param, UInt_t index) const
Check for compatibility between two same-named parameters, which should therefore be identical.
std::vector< Double_t > vectorRes_
Likelihood values returned from the tasks.
Pure abstract base class for defining a parameter containing an R value.
Definition: LauAbsRValue.hh:45
void runSimFit(const TString &fitNtupleFileName, const Bool_t useAsymmErrors=kFALSE, const Bool_t doTwoStageFit=kFALSE)
Run the fit.
TStopwatch cumulTimer_
The total fit timer.
TMonitor * socketMonitor_
Parallel setup monitor.
LauFitNtuple * fitNtuple_
The fit results ntuple.
The coordinator process for simultaneous/combined fits.
virtual Double_t getTotNegLogLikelihood()
Calculate the new value of the negative log likelihood.
LauSimFitCoordinator & operator=(const LauSimFitCoordinator &rhs)
Copy assignment operator (not implemented)
The abstract interface for the objects that control the calculation of the likelihood.
Definition: LauFitObject.hh:47
void checkInitFitParams()
Instruct the tasks to update the initial fit parameter values, if required.
void updateParametersFromTasks()
Update and verify the parameter initial values from all tasks.
const UInt_t nTasks_
The number of tasks.
std::vector< TMessage * > messagesToTasks_
Messages to tasks.
Double_t getLogLikelihoodPenalty()
Calculate the penalty terms to the log likelihood from Gaussian constraints.
File containing declaration of LauFitObject class.
std::vector< TSocket * > socketTasks_
Sockets for each of the tasks.
std::map< UInt_t, TString > parNames_
Reverse map of index in the values vector to parameter names.
void initialise()
Initialise.
Bool_t writeOutResults()
Instruct the tasks to write out the fit results.
std::map< TString, UInt_t > parIndices_
Map of parameter names to index in the values vector.