AliPhysics  2c8507d (2c8507d)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliAnalysisTaskDmesonJetsDetectorResponse.cxx
Go to the documentation of this file.
1 /*************************************************************************
2 * Copyright(c) 1998-2016, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15 
17 
18 // Definitions of class AliAnalysisTaskDmesonJetsDetectorResponse::AliDmesonMatchInfoSummary
19 
23 
24 // Definitions of class AliAnalysisTaskDmesonJetsDetectorResponse::AliD0MatchInfoSummary
25 
29 
32 {
33  fGenerated.Reset();
34  fReconstructed.Reset();
35 }
36 
41 {
42  fReconstructed.Set(reco);
43 }
44 
49 {
50  fGenerated.Set(truth);
51 }
52 
53 // Definitions of class AliAnalysisTaskDmesonJetsDetectorResponse::AliDStarMatchInfoSummary
54 
58 
61 {
62  fGenerated.Reset();
63  fReconstructed.Reset();
64 }
65 
70 {
71  fReconstructed.Set(reco);
72 }
73 
78 {
79  fGenerated.Set(truth);
80 }
81 
82 // Definitions of class AliAnalysisTaskDmesonJetsDetectorResponse::ResponseEngine
83 
87 
90  TObject(),
91  fCandidateType(kD0toKpi),
92  fInhibit(kFALSE),
93  fMaxJetDmesonDistance(0),
94  fName(),
95  fTree(0),
96  fCurrentDmeson(0),
97  fCurrentJetInfoReco(0),
98  fCurrentJetInfoTruth(0),
99  fDataSlotNumber(-1),
100  fRecontructed(0),
101  fGenerated(0)
102 {
103 
104 }
105 
110  TObject(),
111  fCandidateType(type),
112  fInhibit(kFALSE),
113  fMaxJetDmesonDistance(1),
114  fName(),
115  fTree(0),
116  fCurrentDmeson(0),
117  fCurrentJetInfoReco(0),
118  fCurrentJetInfoTruth(0),
119  fDataSlotNumber(-1),
120  fRecontructed(0),
121  fGenerated(0)
122 {
123 
124 }
125 
130  TObject(source),
131  fCandidateType(source.fCandidateType),
132  fInhibit(source.fInhibit),
133  fMaxJetDmesonDistance(1),
134  fName(source.fName),
135  fTree(0),
136  fCurrentDmeson(0),
137  fCurrentJetInfoReco(0),
138  fCurrentJetInfoTruth(0),
139  fDataSlotNumber(source.fDataSlotNumber),
140  fRecontructed(source.fRecontructed),
141  fGenerated(source.fGenerated)
142 {
143 }
144 
149 {
150  new (this) ResponseEngine(source);
151  return *this;
152 }
153 
158 {
159  fInhibit = (fRecontructed == 0 || fGenerated == 0);
160  fName = fRecontructed->GetCandidateName();
161  return !fInhibit;
162 }
163 
166 {
167 
168 }
169 
175 {
176  TString classname;
177  switch (fCandidateType) {
178  case kD0toKpi:
179  case kD0toKpiLikeSign:
180  classname = "AliAnalysisTaskDmesonJetsDetectorResponse::AliD0MatchInfoSummary";
181  fCurrentDmeson = new AliD0MatchInfoSummary();
182  break;
183  case kDstartoKpipi:
184  classname = "AliAnalysisTaskDmesonJetsDetectorResponse::AliDStarMatchInfoSummary";
185  fCurrentDmeson = new AliDStarMatchInfoSummary();
186  break;
187  }
188  TString treeName = TString::Format("%s_%s", taskName, GetName());
189  fTree = new TTree(treeName, treeName);
190  fTree->Branch("DmesonJet", classname, &fCurrentDmeson);
191 
192  fCurrentJetInfoTruth = new AliJetInfoSummary*[fGenerated->GetJetDefinitions().size()];
193  for (Int_t i = 0; i < fGenerated->GetJetDefinitions().size(); i++) {
194  fCurrentJetInfoTruth[i] = new AliJetInfoSummary();
195  TString bname = TString::Format("%s_truth", fGenerated->GetJetDefinitions()[i].GetName());
196  fTree->Branch(bname, "AliAnalysisTaskDmesonJets::AliJetInfoSummary", &fCurrentJetInfoTruth[i]);
197  }
198 
199  fCurrentJetInfoReco = new AliJetInfoSummary*[fRecontructed->GetJetDefinitions().size()];
200  for (Int_t i = 0; i < fRecontructed->GetJetDefinitions().size(); i++) {
201  fCurrentJetInfoReco[i] = new AliJetInfoSummary();
202  TString bname = TString::Format("%s_reco", fRecontructed->GetJetDefinitions()[i].GetName());
203  fTree->Branch(bname, "AliAnalysisTaskDmesonJets::AliJetInfoSummary", &fCurrentJetInfoReco[i]);
204  }
205 
206  return fTree;
207 }
208 
214 {
215  fRecontructed->FillQA(applyKinCuts);
216  fGenerated->FillQA(applyKinCuts);
217 
218  std::map<int, AliDmesonJetInfo>& recoDmesons = fRecontructed->GetDmesons();
219  std::map<int, AliDmesonJetInfo>& truthDmesons = fGenerated->GetDmesons();
220 
221  // Loop over reconstructed D meson and look for their generated counterparts
222  for (auto& dmeson_reco : recoDmesons) {
223  // Reset D meson and jet tree branches
224  fCurrentDmeson->Reset();
225  for (UInt_t ij = 0; ij < fRecontructed->GetJetDefinitions().size(); ij++) {
226  fCurrentJetInfoReco[ij]->Reset();
227  }
228  for (UInt_t ij = 0; ij < fGenerated->GetJetDefinitions().size(); ij++) {
229  fCurrentJetInfoTruth[ij]->Reset();
230  }
231 
232  // Copy the reconstructed D meson information in the tree branch
233  fCurrentDmeson->SetReconstructed(dmeson_reco.second);
234 
235  Int_t accRecJets = 0;
236  Int_t accGenJets = 0;
237 
238  // Copy the reconstructed jet information in the tree branch
239  for (UInt_t ij = 0; ij < fRecontructed->GetJetDefinitions().size(); ij++) {
240  AliJetInfo* jet = dmeson_reco.second.GetJet(fRecontructed->GetJetDefinitions()[ij].GetName());
241  if (!jet) continue;
242  if (applyKinCuts && !fRecontructed->GetJetDefinitions()[ij].IsJetInAcceptance(*jet)) continue;
243  fCurrentJetInfoReco[ij]->Set(dmeson_reco.second, fRecontructed->GetJetDefinitions()[ij].GetName());
244  accRecJets++;
245  }
246 
247  // Look for the generated D meson counterpart
248  if (dmeson_reco.second.fMCLabel >= 0) {
249  std::map<int, AliDmesonJetInfo>::iterator it = truthDmesons.find(dmeson_reco.second.fMCLabel);
250  if (it != truthDmesons.end()) {
251  std::pair<const int, AliDmesonJetInfo>& dmeson_truth = (*it);
252  // Flag the generated D meson as reconstructed
253  dmeson_truth.second.fReconstructed = kTRUE;
254  // Copy the generated D meson information in the tree branch
255  fCurrentDmeson->SetGenerated((*it).second);
256 
257  // Copy the reconstructed jet information in the tree branch
258  for (UInt_t ij = 0; ij < fGenerated->GetJetDefinitions().size(); ij++) {
259  AliJetInfo* jet = dmeson_truth.second.GetJet(fGenerated->GetJetDefinitions()[ij].GetName());
260  if (!jet) continue;
261  if (!applyKinCuts || fGenerated->GetJetDefinitions()[ij].IsJetInAcceptance(*jet)) accGenJets++;
262  fCurrentJetInfoTruth[ij]->Set(dmeson_truth.second, fGenerated->GetJetDefinitions()[ij].GetName());
263  }
264  }
265  }
266 
267  // Fill the tree with the current D meson
268  if (accRecJets > 0 || accGenJets > 0) fTree->Fill();
269  }
270 
271  // Reset jet tree branches
272  for (UInt_t ij = 0; ij < fRecontructed->GetJetDefinitions().size(); ij++) fCurrentJetInfoReco[ij]->Reset();
273  for (UInt_t ij = 0; ij < fGenerated->GetJetDefinitions().size(); ij++) fCurrentJetInfoTruth[ij]->Reset();
274 
275  // Loop over generated D meson that were not reconstructed
276  for (auto& dmeson_truth : truthDmesons) {
277  // Reset D meson tree branch
278  fCurrentDmeson->Reset();
279  // Skip if the generated D meson was reconstructed
280  if (dmeson_truth.second.fReconstructed) continue;
281  // Copy the generated D meson information in the tree branch
282  fCurrentDmeson->SetGenerated(dmeson_truth.second);
283 
284  Int_t accRecJets = 0;
285  Int_t accGenJets = 0;
286 
287  // Copy the reconstructed jet information in the tree branch
288  for (UInt_t ij = 0; ij < fGenerated->GetJetDefinitions().size(); ij++) {
289  fCurrentJetInfoTruth[ij]->Reset();
290  AliJetInfo* jet = dmeson_truth.second.GetJet(fGenerated->GetJetDefinitions()[ij].GetName());
291  if (!jet) continue;
292  if (applyKinCuts && !fGenerated->GetJetDefinitions()[ij].IsJetInAcceptance(*jet)) continue;
293  fCurrentJetInfoTruth[ij]->Set(dmeson_truth.second, fGenerated->GetJetDefinitions()[ij].GetName());
294  accGenJets++;
295  }
296 
297  // Check if a matched reconstructed jet can be found and copy the reconstructed jet information in the tree branch
298  if (findNoDMesonRecoJets) {
299  for (UInt_t ij = 0; ij < fRecontructed->GetJetDefinitions().size(); ij++) {
300  // Reset reconstructed jet tree branch
301  fCurrentJetInfoReco[ij]->Reset();
302  // Look for a matched reconstructed jet
303  AliAnalysisTaskDmesonJets::AnalysisEngine::jet_distance_pair jet = fRecontructed->FindJetMacthedToGeneratedDMeson(dmeson_truth.second, fRecontructed->GetJetDefinitions()[ij], fRecontructed->GetJetDefinitions()[ij].GetRadius() * fMaxJetDmesonDistance, applyKinCuts);
304  if (!jet.first) continue;
305  // Copy reconstructed jet information in the tree branch
306  fCurrentJetInfoReco[ij]->Set(*(jet.first));
307  fCurrentJetInfoReco[ij]->fR = jet.second;
308  accRecJets++;
309  }
310  }
311 
312  // Fill the tree with the current D meson
313  if (accRecJets > 0 || accGenJets > 0) fTree->Fill();
314  }
315 
316  return kTRUE;
317 }
318 
319 // Definitions of class AliAnalysisTaskDmesonJetsDetectorResponse
320 
324 
330 {
332 }
333 
338  AliAnalysisTaskDmesonJets(name, nOutputTrees),
339  fFindRecoJetsForLostDMesons(kFALSE),
340  fResponseEngines()
341 {
343 }
344 
347 {
348  ::Info("UserCreateOutputObjects", "CreateOutputObjects of task %s", GetName());
349 
351 
352  for (auto &param : fAnalysisEngines) {
353  if (param.IsInhibit()) continue;
354  if (param.GetMCMode() != kSignalOnly && param.GetMCMode() != kMCTruth) continue;
355 
356  std::map<ECandidateType_t, ResponseEngine>::iterator it = fResponseEngines.find(param.GetCandidateType());
357 
358  if (it == fResponseEngines.end()) {
359  it = (fResponseEngines.insert(std::pair<const ECandidateType_t, ResponseEngine>(param.GetCandidateType(), ResponseEngine(param.GetCandidateType())))).first;
360  }
361 
362  if (param.GetMCMode() == kMCTruth) {
363  (*it).second.SetGeneratedAnalysisEngine(&param);
364  }
365  else {
366  (*it).second.SetReconstructedAnalysisEngine(&param);
367  }
368  }
369 
370  Int_t treeSlot = 0;
371 
372  for (auto &resp : fResponseEngines) {
373  if (!resp.second.CheckInit()) continue;
374 
375  resp.second.BuildTree(GetName());
376  if (treeSlot < fNOutputTrees) {
377  resp.second.AssignDataSlot(treeSlot+2);
378  treeSlot++;
379  PostDataFromResponseEngine(resp.second);
380  }
381  else {
382  AliError(Form("Number of data output slots %d not sufficient. Tree of response engine %s will not be posted!", fNOutputTrees, resp.second.GetName()));
383  }
384  }
385 }
386 
390 {
392 }
393 
398 {
400 }
401 
406 {
407  for (auto &resp : fResponseEngines) {
408 
409  if (resp.second.IsInhibit()) continue;
410 
411  resp.second.FillTree(fApplyKinematicCuts, fFindRecoJetsForLostDMesons);
412 
413  PostDataFromResponseEngine(resp.second);
414  }
415 
416  return kTRUE;
417 
418 }
419 
424 {
426  AliWarning("This class only provides a tree output.");
427  }
428 
429  // Always set it to kNoOutput: base class does not generate any output (other than QA histograms), all the output comes from the derived class
431 }
432 
439 {
440  if (eng.GetDataSlotNumber() >= 0 && eng.GetTree()) {
441  PostData(eng.GetDataSlotNumber(), eng.GetTree());
442  return eng.GetDataSlotNumber();
443  }
444  else {
445  return -1;
446  }
447 }
std::list< AnalysisEngine > fAnalysisEngines
Array of analysis parameters.
Analysis task for D meson jets.
Lightweight class that encapsulates matching between reconstructed and generated D0 mesons...
Analysis task used to build a detector response for D meson jets.
std::map< ECandidateType_t, ResponseEngine > fResponseEngines
! Response engines
Analysis engine to produce detector response matrix in the D meson jet analysis.
virtual void UserCreateOutputObjects()
Creates the output containers.
virtual void UserCreateOutputObjects()
Creates the output containers.
AliAnalysisTaskDmesonJetsDetectorResponse()
This is the default constructor, used for ROOT I/O purposes.
Lightweight class that encapsulates matching between reconstructed and generated D mesons...
Bool_t fFindRecoJetsForLostDMesons
If switched on, looks for reconstructed jets even when the D meson was lost.
std::pair< AliJetInfo *, Double_t > jet_distance_pair
int Int_t
Definition: External.C:63
unsigned int UInt_t
Definition: External.C:33
Lightweight class that encapsulates matching between reconstructed and generated D* mesons...
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
Bool_t fApplyKinematicCuts
Apply jet kinematic cuts.
EOutputType_t fOutputType
Output type: none, TTree or THnSparse.
bool Bool_t
Definition: External.C:53
Int_t fNOutputTrees
Maximum number of output trees.
void RunAnalysis()
Run the requested analysis for the current event.