3. Analysis in a Task¶
In order to run all the stuff from above as a compiled !FairTask instead of an interpreted ROOT macro some preparations have to be done. But it’s actually easier than you might think! The necessary steps to convert your analysis code to a compilable task are:
Create a class inheriting from FairTask containing the analysis code
Add/modify local files CMakeLists.txt and tutrhoLinkDef.h
Modify CMakeLists.txt in root directory (e.g. trunk/CMakeLists.txt) in order to build the library
Create a macro, which loads the library and adds your analysis task to FairRunAna
Run the macro
Since most of this is not analysis specific but more related to !FairRoot please also take a look to the FairRoot website (http://fairroot.gsi.de/) for further information.
3.1 Create the class PndTutAnaTask¶
In principle the code from our macro tut_ana.C can by copied to the task code almost like it is, but we have to split it up appropriately to the header file PndTutAnaTask.h and the source file PndTutAnaTask.cxx of our class definition. To understand what goes where, let’s first take a look to the structure of a FairTask. The three most important methods of this class are:
Init()
- Initialization of the task. Here the members of the class can be initialized, histograms are create, etc
Exec()
- This is the actual event loop. Here your analysis code is located.
Finish()
- At the end of the analysis, histograms etc are stored.
In order to access variables and objects across the class methods, they have to be declared in the private section of the class definition in the header file PndTutAnaTask.h like it is shown below:
...
class PndTutAnaTask : public FairTask
{
public:
PndTutAnaTask(); // ** Default constructor
~PndTutAnaTask(); // ** Destructor
...
private:
// *** the method from macro
int SelectTruePid(PndAnalysis *ana, RhoCandList &l);
int fEvtCount; // *** event counter
RhoMassParticleSelector *fJpsiMassSel; // *** mass selector for the J/psi
// *** declare some histograms
TH1F *hjpsim_all;
TH1F *hpsim_all;
TH1F *hjpsim_lpid;
...
TLorentzVector fIni; // *** the initial 4-vector
PndAnalysis *fAnalysis; // *** the PndAnalysis object
...
};
...
In particular it is important to declare an object of class PndAnalysis, which works in the task
in the same way as in the macro. After creating the header file, the members like fAnalysis
,
fIni
, the histograms etc have to be setup in the Init()
method in the source file
PndTutAnaTask.cxx:
InitStatus PndTutAnaTask::Init()
{
fAnalysis = new PndAnalysis(); // *** initialize analysis object
fEvtCount = 0; // *** reset the event counter
// *** Mass selector for the jpsi cands
fJpsiMassSel=new RhoMassParticleSelector("jpsi",3.096,1.0);
// *** create the histograms
hjpsim_all = new TH1F("hjpsim_all","J/#psi mass (all)",200,0,4.5);
hpsim_all = new TH1F("hpsim_all","#psi(2S) mass (all)",200,0,5);
...
// *** the lorentz vector of the initial psi(2S)
fIni.SetXYZT(0, 0, 6.231552, 7.240065);
return kSUCCESS;
}
The main part now is to implement the analysis code in the Exec()
method. There we copy
everything from within the event loop from the macro after putting an fAnalysis->GetEventInTask()
(instead of fAnalysis->GetEvent()
used in the macro!) just at the beginning of the method. This
looks like this:
void PndTutAnaTask::Exec(Option_t* opt)
{
// necessary to read the next event
fAnalysis->GetEventInTask(); // instead of GetEvent() in macro
// ---> from now on everything is exactly the same
// ---> as in tut_ana.C, except adaptations of some variable names
// ---> like fIni, fAnalysis,...
// *** RhoCandLists for the analysis
RhoCandList muplus, muminus, piplus, piminus, jpsi, psi2s;
// *** Select with no PID info ('All'); type and mass are set
fAnalysis->FillList(muplus, "MuonAllPlus");
fAnalysis->FillList(muminus, "MuonAllMinus");
fAnalysis->FillList(piplus, "PionAllPlus");
fAnalysis->FillList(piminus, "PionAllMinus");
// *** combinatorics for J/psi -> mu+ mu-
jpsi.Combine(muplus, muminus);
...
}
Finally in the Finish()
method, all the histograms are stored in the output file. Since the IO is
handled by FairRootManager, we don’t have to care about opening files etc. The corresponding code
looks like this:
void PndTutAnaTask::Finish()
{
// *** Store all the histograms
hjpsim_all->Write();
hpsim_all->Write();
hjpsim_lpid->Write();
...
}
3.2 Compile the class and build a library¶
The code of our new task now is ready, we just need to tell PandaROOT to compile it and build a corresponding library which we can use in a macro. Therefore we basically need two local files: CMakeLists.txt and tutrhoLinkDef.h. The first one is a
configuration file for cmake
containing all dependencies of your class (i.e. from which other !PandaROOT directories it needs to include stuff), the name and location of source and header files and the name of the library to be built. The most important part here is to add your source and header files (in this example the library will be named libtutrhotask.so) in CMakeLists.txt:
set(tutrhotask_SRCS
...
PndTutAnaTask.cxx
)
set(tutrhotask_HEADERS
...
PndTutAnaTask.h
In addition you need to add a line in the tutrhoLinkDef.h file like this:
#pragma link C++ class PndTutAnaTask+;
Last but not least, you have to modify the CMakeLists.txt in your $VMCWORKDIR
(e.g. trunk/) by adding a line:
add_subdirectory(tutorials/rho)
somewhere close to the place, where all the other add_subdirectory
commands are.
Now you are able to compile your code and build the library by just executing make
in your build directory.
3.3 Loading library and running the macro¶
For running the analysis we still need a macro, but in contrast to tut_ana.C it will look fairly simple and is completely independent of the actual analysis code. It starts like every standard PandaROOT macro, where you initialize some things and define your input and output files. You then just need to add your task to the FairRunAna and start the event loop. The full macro looks like this:
void tut_ana_task(int nevts = 0, TString prefix = "signal")
{
TString OutFile=prefix+"_ana_task.root";
// *** the files coming from the simulation
TString inPidFile = prefix+"_pid.root"; // PndPidCandidates and McTruth
TString inParFile = prefix+"_par.root";
// *** PID table with selection thresholds; can be modified by the user
TString pidParFile = TString(gSystem->Getenv("VMCWORKDIR"))+"/macro/params/all.par";
// *** initialization
FairLogger::GetLogger()->SetLogToFile(kFALSE);
FairRunAna* fRun = new FairRunAna();
FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
fRun->SetSource(new FairFileSource(inPidFile));
// *** setup parameter database
FairParRootFileIo* parIO = new FairParRootFileIo();
parIO->open(inParFile);
FairParAsciiFileIo* parIOPid = new FairParAsciiFileIo();
parIOPid->open(pidParFile.Data(),"in");
rtdb->setFirstInput(parIO);
rtdb->setSecondInput(parIOPid);
rtdb->setOutput(parIO);
fRun->SetOutputFile(OutFile);
// *** HERE OUR TASK GOES!
PndTutAnaTask *anaTask = new PndTutAnaTask();
fRun->AddTask(anaTask);
// *** and run analysis
fRun->Init();
fRun->Run(0,nevts);
}
Now you can run your macro with root -l -b -q tut_ana_task.C
, and the analysis should just run like with tut_ana.C.