FairRoot/PandaRoot
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PndMvdRadDamTask Class Reference

#include <PndMvdRadDamTask.h>

Inheritance diagram for PndMvdRadDamTask:

Public Member Functions

 PndMvdRadDamTask ()
 
 ~PndMvdRadDamTask ()
 
 PndMvdRadDamTask (const PndMvdRadDamTask &)=delete
 
PndMvdRadDamTaskoperator= (const PndMvdRadDamTask &)=delete
 
virtual void SetParContainers ()
 
virtual InitStatus Init ()
 
virtual InitStatus ReInit ()
 
virtual void Exec (Option_t *opt)
 
virtual void Finish ()
 
void SetPersistance (Bool_t p=kTRUE)
 
Bool_t GetPersistance ()
 

Private Member Functions

void InitWeightLists ()
 
 ClassDef (PndMvdRadDamTask, 1)
 

Private Attributes

Bool_t fPersistance
 
TClonesArray * fMCTracks
 
TClonesArray * fMCHits
 
TClonesArray * fRadDamHits
 
PndMvdRadDamListfElectronList
 
PndMvdRadDamListfProtonList
 
PndMvdRadDamListfNeutronList
 
PndMvdRadDamListfPionList
 
PndGeoHandlingfGeoH
 
std::map< std::string, TH2 * > fMapDetHistos
 
TH1D * fRadDamHisto
 
std::map< Int_t,
PndMvdRadDamList * > 
fWeightListsMap
 

Detailed Description

Definition at line 21 of file PndMvdRadDamTask.h.

Constructor & Destructor Documentation

PndMvdRadDamTask::PndMvdRadDamTask ( )

Definition at line 19 of file PndMvdRadDamTask.cxx.

19  :
20  fPersistance(kTRUE),
21  fMCTracks(NULL),
22  fMCHits(NULL),
23  fRadDamHits(NULL),
24  fElectronList(NULL),
25  fProtonList(NULL),
26  fNeutronList(NULL),
27  fPionList(NULL),
29  fMapDetHistos(),
30  fRadDamHisto(NULL),
32 {
33 }
TClonesArray * fMCHits
TClonesArray * fMCTracks
std::map< Int_t, PndMvdRadDamList * > fWeightListsMap
PndGeoHandling * fGeoH
PndMvdRadDamList * fNeutronList
TClonesArray * fRadDamHits
static PndGeoHandling * Instance()
PndMvdRadDamList * fProtonList
PndMvdRadDamList * fPionList
PndMvdRadDamList * fElectronList
std::map< std::string, TH2 * > fMapDetHistos
PndMvdRadDamTask::~PndMvdRadDamTask ( )

Definition at line 35 of file PndMvdRadDamTask.cxx.

References fMapDetHistos, and fRadDamHisto.

36 {
37  for (std::map<std::string, TH2*>::iterator it = fMapDetHistos.begin(); it != fMapDetHistos.end(); it++)
38  delete(it->second);
39  fMapDetHistos.clear();
40  delete (fRadDamHisto);
41 }
std::map< std::string, TH2 * > fMapDetHistos
PndMvdRadDamTask::PndMvdRadDamTask ( const PndMvdRadDamTask )
delete

Member Function Documentation

PndMvdRadDamTask::ClassDef ( PndMvdRadDamTask  ,
 
)
private
void PndMvdRadDamTask::Exec ( Option_t *  opt)
virtual

Virtual method Exec

Definition at line 115 of file PndMvdRadDamTask.cxx.

References detname, Double_t, fGeoH, fMapDetHistos, fMCHits, fMCTracks, fRadDamHisto, fRadDamHits, fWeightListsMap, PndGeoHandling::GetPath(), PndMCTrack::GetPdgCode(), PndSdsMCPoint::GetPosition(), PndGeoHandling::GetSensorDimensionsShortId(), PndSdsMCPoint::GetSensorID(), i, PndGeoHandling::MasterToLocalShortId(), mom, pid(), and TString.

116 {
117  // Reset output array
118  if ( ! fRadDamHits )
119  Fatal("Exec", "No RadDamArray");
120  fRadDamHits->Delete();
121 
122  PndSdsMCPoint* mcPoint;
123  PndMCTrack* mcTrack;
124 
125  for (int i = 0; i < fMCHits->GetEntriesFast(); i++){
126  mcPoint = (PndSdsMCPoint*)fMCHits->At(i);
127  mcTrack = (PndMCTrack*)(fMCTracks->At(mcPoint->GetTrackID()));
128  TVector3 mom;
129  mcPoint->Momentum(mom);
130  Double_t Ekin = mom.Mag();
131  Int_t pid = mcTrack->GetPdgCode();
132  Double_t weight = 0;
133  if (fWeightListsMap[pid] != NULL){
134  weight = fWeightListsMap[pid]->GetWeight(Ekin);
135  }
136  //std::cout << "WeightCalc: pid: " << pid << " Energy: " << Ekin << " Weight: " << weight << std::endl;
137  new ((*fRadDamHits)[i]) PndMvdRadDamHit(mcPoint->GetTrackID(), i, mcPoint->GetSensorID(), pid, Ekin,
138  mcPoint->GetPosition(), mom, weight);
139  TString detname = fGeoH->GetPath(mcPoint->GetSensorID());
140  if (fMapDetHistos[detname.Data()] == 0){
141  PndStringSeparator svec(detname.Data(),"/");
142  TVector3 sensDim = fGeoH->GetSensorDimensionsShortId(mcPoint->GetSensorID());
143  fMapDetHistos[detname.Data()] = new TH2D(svec.Replace("/","o").c_str(),
144  fGeoH->GetPath(mcPoint->GetSensorID()),
145  (Int_t)(2*sensDim.X()*10),-sensDim.X(),sensDim.X(), // point resolution mm^2
146  (Int_t)(2*sensDim.Y()*10),-sensDim.Y(),sensDim.Y());
147  }
148  TVector3 localHit = fGeoH->MasterToLocalShortId(mcPoint->GetPosition(), mcPoint->GetSensorID());
149  ((TH2D*)(fMapDetHistos[detname.Data()]))->Fill(localHit.X(), localHit.Y(), weight);
150  fRadDamHisto->Fill(weight);
151  }
152 }
TClonesArray * fMCHits
Int_t i
Definition: run_full.C:25
TClonesArray * fMCTracks
std::map< Int_t, PndMvdRadDamList * > fWeightListsMap
Int_t GetPdgCode() const
Definition: PndMCTrack.h:73
PndGeoHandling * fGeoH
int pid()
TString detname
Definition: anasim.C:61
Double_t mom
Definition: plot_dirc.C:14
Int_t GetSensorID() const
Definition: PndSdsMCPoint.h:89
TString GetPath(Int_t shortID)
for a given shortID the path is returned
TVector3 GetSensorDimensionsShortId(Int_t shortId)
Double_t
TClonesArray * fRadDamHits
TVector3 MasterToLocalShortId(const TVector3 &master, const Int_t &shortId)
TVector3 GetPosition() const
Definition: PndSdsMCPoint.h:90
std::map< std::string, TH2 * > fMapDetHistos
void PndMvdRadDamTask::Finish ( )
virtual

Definition at line 154 of file PndMvdRadDamTask.cxx.

References fMapDetHistos, and fRadDamHisto.

155 {
156  for (std::map<std::string, TH2*>::iterator it = fMapDetHistos.begin(); it != fMapDetHistos.end(); it++)
157  it->second->Write();
158  fRadDamHisto->Write();
159 }
std::map< std::string, TH2 * > fMapDetHistos
Bool_t PndMvdRadDamTask::GetPersistance ( )
inline

Definition at line 39 of file PndMvdRadDamTask.h.

References fPersistance.

39 {return fPersistance;};
InitStatus PndMvdRadDamTask::Init ( )
virtual

Definition at line 54 of file PndMvdRadDamTask.cxx.

References fMCHits, fMCTracks, fPersistance, fRadDamHisto, fRadDamHits, and InitWeightLists().

55 {
56  FairRootManager* ioman = FairRootManager::Instance();
57 
58  if (!ioman)
59  {
60  std::cout << "-E- PndMvdRadDamTask::Init: "
61  << "RootManager not instantiated!" << std::endl;
62  return kFATAL;
63  }
64 
65  fMCTracks = (TClonesArray*) ioman->GetObject("MCTrack");
66  if (!fMCTracks)
67  {
68  std::cout << "-W- PndMvdRadDamTask::Init: " << "No MCTrack array!"
69  << std::endl;
70  return kERROR;
71  }
72 
73  fMCHits = (TClonesArray*) ioman->GetObject("MVDPoint");
74  if (!fMCHits)
75  {
76  std::cout << "-W- PndMvdRadDamTask::Init: " << "No MVDPoint array!"
77  << std::endl;
78  return kERROR;
79  }
80 
81  // Create and register output array
82  fRadDamHits = new TClonesArray("PndMvdRadDamHit");
83  ioman->Register("MVDRadDamHit", "MVD", fRadDamHits, fPersistance);
84 
86  fRadDamHisto = new TH1D("radDamH","Weight Factors", 1000,0,100);
87 
88  std::cout << "-I- PndMvdRadDamTask: Initialization successful" << std::endl;
89 
90  return kSUCCESS;
91 }
TClonesArray * fMCHits
TClonesArray * fMCTracks
TClonesArray * fRadDamHits
void PndMvdRadDamTask::InitWeightLists ( )
private

Definition at line 93 of file PndMvdRadDamTask.cxx.

References fElectronList, fNeutronList, fPionList, fProtonList, and fWeightListsMap.

Referenced by Init().

94 {
95  fElectronList = new PndMvdRadDamList("$VMCWORKDIR/pandaroot/mvd/MvdTools/MvdRadDamage/electronsWeight.root");
96  fProtonList = new PndMvdRadDamList("$VMCWORKDIR/pandaroot/mvd/MvdTools/MvdRadDamage/protonsWeight.root");
97  fNeutronList = new PndMvdRadDamList("$VMCWORKDIR/pandaroot/mvd/MvdTools/MvdRadDamage/neutronsWeight.root");
98  fPionList = new PndMvdRadDamList("$VMCWORKDIR/pandaroot/mvd/MvdTools/MvdRadDamage/pionsWeight.root");
99 
100  fWeightListsMap[11] = fElectronList; //e-
101  fWeightListsMap[-11] = fElectronList; //e+
102  fWeightListsMap[2212] = fProtonList; //p+
103  fWeightListsMap[-2212] = fProtonList; //p-
104  fWeightListsMap[2112] = fNeutronList; //n
105  fWeightListsMap[-2112] = fNeutronList; //n_bar
106  fWeightListsMap[111] = fPionList; //pi0
107  fWeightListsMap[211] = fPionList; //pi+
108  fWeightListsMap[-211] = fPionList; //pi-
109  fWeightListsMap[1000010030] = fProtonList; //Tritium
110  fWeightListsMap[1000010020] = fProtonList; //Deuterium
111 
112 
113 }
std::map< Int_t, PndMvdRadDamList * > fWeightListsMap
PndMvdRadDamList * fNeutronList
PndMvdRadDamList * fProtonList
PndMvdRadDamList * fPionList
PndMvdRadDamList * fElectronList
PndMvdRadDamTask& PndMvdRadDamTask::operator= ( const PndMvdRadDamTask )
delete
InitStatus PndMvdRadDamTask::ReInit ( )
virtual

Definition at line 48 of file PndMvdRadDamTask.cxx.

References SetParContainers().

49 {
51  return kSUCCESS;
52 }
virtual void SetParContainers()
void PndMvdRadDamTask::SetParContainers ( )
virtual

Definition at line 44 of file PndMvdRadDamTask.cxx.

Referenced by ReInit().

45 {
46 }
void PndMvdRadDamTask::SetPersistance ( Bool_t  p = kTRUE)
inline

Definition at line 38 of file PndMvdRadDamTask.h.

References fPersistance, and p.

38 {fPersistance=p;};
Double_t p
Definition: anasim.C:58

Member Data Documentation

PndMvdRadDamList* PndMvdRadDamTask::fElectronList
private

Definition at line 47 of file PndMvdRadDamTask.h.

Referenced by InitWeightLists().

PndGeoHandling* PndMvdRadDamTask::fGeoH
private

Definition at line 52 of file PndMvdRadDamTask.h.

Referenced by Exec().

std::map<std::string, TH2*> PndMvdRadDamTask::fMapDetHistos
private

Definition at line 54 of file PndMvdRadDamTask.h.

Referenced by Exec(), Finish(), and ~PndMvdRadDamTask().

TClonesArray* PndMvdRadDamTask::fMCHits
private

Definition at line 44 of file PndMvdRadDamTask.h.

Referenced by Exec(), and Init().

TClonesArray* PndMvdRadDamTask::fMCTracks
private

Definition at line 43 of file PndMvdRadDamTask.h.

Referenced by Exec(), and Init().

PndMvdRadDamList* PndMvdRadDamTask::fNeutronList
private

Definition at line 49 of file PndMvdRadDamTask.h.

Referenced by InitWeightLists().

Bool_t PndMvdRadDamTask::fPersistance
private

Definition at line 39 of file PndMvdRadDamTask.h.

Referenced by GetPersistance(), Init(), and SetPersistance().

PndMvdRadDamList* PndMvdRadDamTask::fPionList
private

Definition at line 50 of file PndMvdRadDamTask.h.

Referenced by InitWeightLists().

PndMvdRadDamList* PndMvdRadDamTask::fProtonList
private

Definition at line 48 of file PndMvdRadDamTask.h.

Referenced by InitWeightLists().

TH1D* PndMvdRadDamTask::fRadDamHisto
private

Definition at line 55 of file PndMvdRadDamTask.h.

Referenced by Exec(), Finish(), Init(), and ~PndMvdRadDamTask().

TClonesArray* PndMvdRadDamTask::fRadDamHits
private

Definition at line 45 of file PndMvdRadDamTask.h.

Referenced by Exec(), and Init().

std::map<Int_t, PndMvdRadDamList*> PndMvdRadDamTask::fWeightListsMap
private

Definition at line 57 of file PndMvdRadDamTask.h.

Referenced by Exec(), and InitWeightLists().


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