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

splits clusters based on local maxima in the Phi direction for use with Bremstrahlung correction. More...

#include <PndEmcPhiBumpSplitter.h>

Inheritance diagram for PndEmcPhiBumpSplitter:
PndPersistencyTask

Public Member Functions

 PndEmcPhiBumpSplitter (Int_t verbose=0)
 
virtual ~PndEmcPhiBumpSplitter ()
 
virtual InitStatus Init ()
 Init Task. More...
 
virtual void Exec (Option_t *opt)
 Runs the task. More...
 
virtual void FinishTask ()
 Called at end of task. More...
 
void SetStorageOfData (Bool_t p=kTRUE)
 
PndEmcBumpAddPhiBump ()
 Adds a new PndEmcBump to fPhiBumpArray and returns it. More...
 
void SetPersistency (Bool_t val=kTRUE)
 
Bool_t GetPersistency ()
 

Protected Member Functions

virtual void SetParContainers ()
 

Private Member Functions

 PndEmcPhiBumpSplitter (const PndEmcPhiBumpSplitter &)
 
PndEmcPhiBumpSplitteroperator= (const PndEmcPhiBumpSplitter &)
 
 ClassDef (PndEmcPhiBumpSplitter, 2)
 

Private Attributes

TClonesArray * fDigiArray
 
TClonesArray * fClusterArray
 
TClonesArray * fPhiBumpArray
 
PndEmcGeoParfGeoPar
 
PndEmcDigiParfDigiPar
 
PndEmcRecoParfRecoPar
 
std::vector< Double_tfClusterPosParam
 

Detailed Description

splits clusters based on local maxima in the Phi direction for use with Bremstrahlung correction.

Definition at line 61 of file PndEmcPhiBumpSplitter.h.

Constructor & Destructor Documentation

PndEmcPhiBumpSplitter::PndEmcPhiBumpSplitter ( Int_t  verbose = 0)

Definition at line 75 of file PndEmcPhiBumpSplitter.cxx.

References fClusterPosParam, and PndPersistencyTask::SetPersistency().

75  :PndPersistencyTask("PndEmcPhiBumpSplitter", verbose),
77 
78 {
79  fClusterPosParam.clear();
80  SetPersistency(kTRUE);
81 }
#define verbose
void SetPersistency(Bool_t val=kTRUE)
std::vector< Double_t > fClusterPosParam
parameter set of Emc digitisation
Definition: PndEmcDigiPar.h:12
Parameter set for Emc Reco.
Definition: PndEmcRecoPar.h:12
PndEmcPhiBumpSplitter::~PndEmcPhiBumpSplitter ( )
virtual

Definition at line 87 of file PndEmcPhiBumpSplitter.cxx.

88 {
89 }
PndEmcPhiBumpSplitter::PndEmcPhiBumpSplitter ( const PndEmcPhiBumpSplitter )
private

Member Function Documentation

PndEmcBump * PndEmcPhiBumpSplitter::AddPhiBump ( )

Adds a new PndEmcBump to fPhiBumpArray and returns it.

Returns
PndEmcBump*

Definition at line 283 of file PndEmcPhiBumpSplitter.cxx.

References fPhiBumpArray.

Referenced by Exec().

284 {
285  TClonesArray& clref = *fPhiBumpArray;
286  Int_t size = clref.GetEntriesFast();
287  return new(clref[size]) PndEmcBump();
288 }
represents a reconstructed (splitted) emc cluster
Definition: PndEmcBump.h:34
PndEmcPhiBumpSplitter::ClassDef ( PndEmcPhiBumpSplitter  ,
 
)
private
void PndEmcPhiBumpSplitter::Exec ( Option_t *  opt)
virtual

Runs the task.

Parameters
optunused
Returns
void

Definition at line 149 of file PndEmcPhiBumpSplitter.cxx.

References AddPhiBump(), PndEmcCluster::DigiList(), Double_t, fabs(), fClusterArray, fDigiArray, fPhiBumpArray, PndEmcDigi::GetEnergy(), PndEmcDigi::GetPhi(), i, PndEmcMapper::Instance(), PndEmcBump::MadeFrom(), p, PndEmcCluster::position(), rotate(), PndEmcCluster::SetEnergy(), and PndEmcCluster::SetPosition().

150 {
151  PndEmcMapper::Instance(); //PndEmcMapper *fEmcMap= //[R.K.03/2017] unused variable
152 
153  // Reset output array
154  if ( ! fPhiBumpArray ) Fatal("Exec", "No Phi-Bump Array");
155  fPhiBumpArray->Delete();
156 
157  // loop on each cluster. For each cluster there can be any number of phi_bumps (atleast one)
158  int nClusters = fClusterArray->GetEntriesFast();
159  for (Int_t iCluster = 0; iCluster < nClusters; iCluster++){
160 
161  PndEmcCluster* theCluster = (PndEmcCluster*) fClusterArray->At(iCluster);
162 
163  const Int_t TotNumOfHitPhi = 160;
164 
165  std::vector<double> phi_bump(TotNumOfHitPhi,0);
166 
167  Int_t digiSize = theCluster->DigiList().size();
168  std::vector<Int_t> digiList = theCluster->DigiList();
169  for (Int_t i_digi = 0; i_digi<digiSize; ++i_digi)
170  {
171  PndEmcDigi *emcDigi = (PndEmcDigi *) fDigiArray->At(digiList[i_digi]);
172  Double_t emcDigiPhi = emcDigi->GetPhi()*TMath::RadToDeg();
173  Double_t emcDigiEnergy = emcDigi->GetEnergy();
174  if ( fabs(emcDigiPhi)<=180 ) {
175  Int_t iEmcDigiPhi = int( (emcDigiPhi + 180.) * TotNumOfHitPhi / 360. );
176  phi_bump.at(iEmcDigiPhi) += emcDigiEnergy;
177  }
178  }
179 
180  std::vector<double> vPhiList;
181  std::vector<double> vDepoEnergyList;
182  std::vector<int> vGapSizeList;
183 
184  int i_phi_prev = 0;
185  for (Int_t i_phi=0; i_phi < TotNumOfHitPhi; ++i_phi) {
186  Double_t BinValue = phi_bump.at(i_phi);
187  if (BinValue != 0)
188  {
189  vPhiList.push_back(-180. + ((0.5+i_phi)*360./TotNumOfHitPhi));
190  vDepoEnergyList.push_back(BinValue);
191  vGapSizeList.push_back(i_phi - i_phi_prev );
192  i_phi_prev = i_phi;
193  }
194  }
195 
196  // Find start bin number for the phi projection of cluster.
197  Int_t StartIndex = 0;
198  for (size_t i = 0;i < vGapSizeList.size();i++)
199  {
200  if (vGapSizeList.at(i) > vGapSizeList.at(StartIndex)) StartIndex = i;
201  }
202 
203  // Rotate to the begining of the cluster in case cluster folds over to index=0
204  std::rotate(vDepoEnergyList.begin(),vDepoEnergyList.begin()+StartIndex,vDepoEnergyList.end());
205  vDepoEnergyList.push_back(0);
206  vDepoEnergyList.push_back(0);
207  std::rotate(vDepoEnergyList.begin(),vDepoEnergyList.begin()+(vDepoEnergyList.size()-1),vDepoEnergyList.end());
208  // Do the same for the phi positions
209  std::rotate(vPhiList.begin(),vPhiList.begin()+StartIndex,vPhiList.end());
210  vPhiList.push_back(0);
211  vPhiList.push_back(0);
212  std::rotate(vPhiList.begin(),vPhiList.begin()+(vPhiList.size()-1),vPhiList.end());
213 
214  // Loop through deposited energy vector and identify "valley" type bins => -_- and calculate wieghts to split energy
215  std::vector<int> ValleyType;
216  std::vector<double> Weight;
217  double _Weight = 0;
218  Weight.push_back(0);
219  ValleyType.push_back(0);
220  for (size_t n_sel = 1; n_sel < vDepoEnergyList.size()-1; n_sel++)
221  {
222  if(vDepoEnergyList.at(n_sel-1) > vDepoEnergyList.at(n_sel) &&
223  vDepoEnergyList.at(n_sel) < vDepoEnergyList.at(n_sel+1) ) {
224  ValleyType.push_back(1);
225  Weight.push_back(_Weight);
226  _Weight = 0;
227  } else {
228  _Weight += vDepoEnergyList.at(n_sel);
229  ValleyType.push_back(0);
230  }
231  }
232  Weight.push_back(_Weight);
233  Weight.push_back(0);
234 
235  std::vector<double> enePhiBump, phiPhiBump;
236  int ValleyIndex = 0;
237  int iWeight = 0;
238  for (size_t n_sel = 1; n_sel < vDepoEnergyList.size()-1; n_sel++)
239  {
240  if (ValleyType.at(n_sel) == 1 || n_sel == vDepoEnergyList.size()-2)
241  {
242  iWeight++;
243 
244  const double _eneRightEdge = vDepoEnergyList.at(ValleyIndex)*(Weight.at(iWeight)/(Weight.at(iWeight)+Weight.at(iWeight-1)));
245  double _enePhiBump = _eneRightEdge;
246  double _phiPhiBump = vPhiList.at(ValleyIndex)*_eneRightEdge;
247  for(size_t p = ValleyIndex+1;p < n_sel;p++) {
248  _enePhiBump += vDepoEnergyList.at(p);
249  _phiPhiBump += vPhiList.at(p)*vDepoEnergyList.at(p);
250  }
251  const double _eneLeftEdge = vDepoEnergyList.at(n_sel)*(Weight.at(iWeight)/(Weight.at(iWeight)+Weight.at(iWeight+1) ));
252  _enePhiBump += _eneLeftEdge;
253  _phiPhiBump += vPhiList.at(n_sel)*_eneLeftEdge;
254  _phiPhiBump /= _enePhiBump;
255  enePhiBump.push_back(_enePhiBump);
256  phiPhiBump.push_back(_phiPhiBump);
257  ValleyIndex = n_sel;
258  }
259  }
260 
261  TVector3 posClust = theCluster->position();
262  for (size_t i_phibump=0; i_phibump<enePhiBump.size(); ++i_phibump) {
263  PndEmcBump* theNewPhiBump = AddPhiBump();
264  theNewPhiBump->MadeFrom(iCluster);
265  theNewPhiBump->SetLink(FairLink("EmcCluster", iCluster));
266  theNewPhiBump->SetEnergy(enePhiBump.at(i_phibump));
267  theNewPhiBump->SetTimeStamp(theCluster->GetTimeStamp());
268  theNewPhiBump->SetTimeStampError(theCluster->GetTimeStampError());
269  TVector3 posPhiBump;
270  posPhiBump.SetMagThetaPhi(posClust.Mag(),posClust.Theta(),phiPhiBump.at(i_phibump)*TMath::DegToRad());
271  theNewPhiBump->SetPosition(posPhiBump);
272  }
273 
274  }
275 
276 }
void SetEnergy(Double_t en)
virtual void MadeFrom(Int_t clusterIndex)
Definition: PndEmcBump.cxx:76
Double_t p
Definition: anasim.C:58
virtual Double_t GetEnergy() const
Definition: PndEmcDigi.cxx:296
represents the reconstructed hit of one emc crystal
Definition: PndEmcDigi.h:40
Int_t i
Definition: run_full.C:25
const std::vector< Int_t > & DigiList() const
Definition: PndEmcCluster.h:40
TVector3 position() const
Double_t
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
a cluster (group of neighboring crystals) of hit emc crystals
Definition: PndEmcCluster.h:29
TVector3 rotate(TVector3 vec, TGeoRotation *rotma)
Double_t GetPhi() const
Definition: PndEmcDigi.h:102
void SetPosition(TVector3 pos)
static PndEmcMapper * Instance()
represents a reconstructed (splitted) emc cluster
Definition: PndEmcBump.h:34
PndEmcBump * AddPhiBump()
Adds a new PndEmcBump to fPhiBumpArray and returns it.
void PndEmcPhiBumpSplitter::FinishTask ( )
virtual

Called at end of task.

Returns
void

Definition at line 295 of file PndEmcPhiBumpSplitter.cxx.

296 {
297  cout<<"================================================="<<endl;
298  cout<<"PndEmcPhiBumpSplitter::FinishTask"<<endl;
299  cout<<"================================================="<<endl;
300 }
Bool_t PndPersistencyTask::GetPersistency ( )
inlineinherited

Definition at line 32 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency.

Referenced by PndLmdPixelHitProducerFast::GetPersistance(), PndMdtDigitization::Init(), PndMdtHitProducerIdeal::Init(), PndMdtClusterTask::Init(), PndFtsHitProducerRealFast::Init(), PndDiscTaskReconstruction::Init(), PndRichHitProducer::Init(), PndSttHitProducerRealFast::Init(), PndSttHelixHitProducer::Init(), PndDiscTaskPID::Init(), PndIdealTrackFinder::Init(), PndSttMvdGemTracking::Init(), PndMdtTrkProducer::Init(), PndFtsHitProducerRealFull::Init(), PndLmdPixelClusterTask::Init(), PndSttHitProducerRealFull::Init(), PndLmdStripClusterTask::Init(), PndEmcApdHitProducer::Init(), PndMissingPzCleanerTask::Init(), PndEmcMakeRecoHit::Init(), PndEmcMakeClusterOnline::Init(), PndTrackSmearTask::Init(), PndEmcFWEndcapTimebasedWaveforms::Init(), PndSttHitProducerIdeal::Init(), PndEmcFWEndcapDigi::Init(), PndFtsHitProducerIdeal::Init(), PndEmcMakeCluster::Init(), PndMdtPointsToWaveform::Init(), PndDiscTaskDigitization::Init(), PndEmcMakeDigi::Init(), PndSdsTimeWalkCorrTask::Init(), PndLmdPixelHitProducerFast::Init(), PndDrcHitFinder::Init(), PndRichHitFinder::Init(), PndEmcMakeCorr::Init(), PndFtofHitProducerIdeal::Init(), PndEmcHitsToWaveform::Init(), PndSciTDigiTask::Init(), PndDrcHitProducerIdeal::Init(), PndSdsHitProducerIdeal::Init(), PndSciTHitProducerIdeal::Init(), PndRecoMultiKalmanTask2::Init(), PndEmcHitProducer::Init(), PndDrcHitProducerReal::Init(), PndDskFLGHitProducerIdeal::Init(), PndEmcTmpWaveformToDigi::Init(), PndDrcDigiTask::Init(), PndEmcWaveformToDigi::Init(), PndSttMatchTracks::Init(), PndEmcWaveformToCalibratedDigi::Init(), PndTrkTracking2::Init(), PndSttFindTracks::Init(), PndEmcMultiWaveformToCalibratedDigi::Init(), PndRecoKalmanTask2::Init(), PndDrcTimeDigiTask::Init(), PndEmcExpClusterSplitter::Init(), PndFtsHoughTrackerTask::Init(), PndSdsNoiseProducer::Init(), Init(), PndSdsIdealRecoTask::Init(), PndSdsHybridHitProducer::Init(), PndRecoMultiKalmanTask::Init(), PndSdsIdealClusterTask::Init(), PndRecoKalmanTask::Init(), PndSdsStripHitProducerDif::Init(), PndGemDigitize::Init(), PndSdsStripHitProducer::Init(), PndGemFindHits::Init(), PndSdsPixelClusterTask::Init(), PndSdsStripClusterTask::Init(), PndMvdGemTrackFinderOnHits::Init(), PndBarrelTrackFinder::Init(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcMakeBump::PndEmcMakeBump(), PndUnassignedHitsTask::RegisterBranches(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndEmcMakeBump::SetStorageOfData(), and PndEmcFullDigiTask::StoreDigi().

32 { return fPersistency; }
InitStatus PndEmcPhiBumpSplitter::Init ( )
virtual

Init Task.

Returns
InitStatus
Return values
kSUCCESSsuccess

Definition at line 97 of file PndEmcPhiBumpSplitter.cxx.

References fClusterArray, fClusterPosParam, fDigiArray, fGeoPar, fPhiBumpArray, fRecoPar, PndEmcRecoPar::GetEmcClusterPosMethod(), PndEmcRecoPar::GetOffsetParmA(), PndEmcRecoPar::GetOffsetParmB(), PndEmcRecoPar::GetOffsetParmC(), PndPersistencyTask::GetPersistency(), PndEmcGeoPar::InitEmcMapper(), and PndEmcStructure::Instance().

98 {
99  // Get RootManager
100  FairRootManager* ioman = FairRootManager::Instance();
101  if ( ! ioman ){
102  cout << "-E- PndEmcPhiBumpSplitter::Init: "
103  << "RootManager not instantiated!" << endl;
104  return kFATAL;
105  }
106 
107  // Geometry loading
110 
111  // Get input array
112  fDigiArray = dynamic_cast<TClonesArray *> (ioman->GetObject("EmcDigi"));
113  if ( ! fDigiArray ) {
114  cout << "-W- PndEmcPhiBumpSplitter::Init: "
115  << "No PndEmcDigi array!" << endl;
116  return kERROR;
117  }
118 
119  fClusterArray = dynamic_cast<TClonesArray *> (ioman->GetObject("EmcCluster"));
120  if ( ! fClusterArray ) {
121  cout << "-W- PndEmcPhiBumpSplitter::Init: "
122  << "No PndEmcCluster array!" << endl;
123  return kERROR;
124  }
125 
126  // Set minimum SharedDigi energy to 20keV.
127  if (!strcmp(fRecoPar->GetEmcClusterPosMethod(),"lilo")) {
128  cout<<"Lilo cluster position method"<<endl;
132  }
133 
134  // Create and register output array
135  fPhiBumpArray = new TClonesArray("PndEmcBump");
136  ioman->Register("EmcPhiBump","Emc",fPhiBumpArray,GetPersistency());
137 
138  cout << "-I- PndEmcPhiBumpSplitter: Intialization successfull" << endl;
139 
140  return kSUCCESS;
141 }
std::vector< Double_t > fClusterPosParam
Double_t GetOffsetParmB()
Definition: PndEmcRecoPar.h:22
Double_t GetOffsetParmC()
Definition: PndEmcRecoPar.h:23
void InitEmcMapper()
Double_t GetOffsetParmA()
Definition: PndEmcRecoPar.h:21
Text_t * GetEmcClusterPosMethod()
Definition: PndEmcRecoPar.h:20
static PndEmcStructure * Instance()
PndEmcPhiBumpSplitter& PndEmcPhiBumpSplitter::operator= ( const PndEmcPhiBumpSplitter )
private
void PndEmcPhiBumpSplitter::SetParContainers ( )
protectedvirtual

Get parameter containers

Definition at line 302 of file PndEmcPhiBumpSplitter.cxx.

References fDigiPar, fGeoPar, fRecoPar, and run.

302  {
303 
304  // Get run and runtime database
305  FairRun* run = FairRun::Instance();
306  if ( ! run ) Fatal("SetParContainers", "No analysis run");
307 
308  FairRuntimeDb* db = run->GetRuntimeDb();
309  if ( ! db ) Fatal("SetParContainers", "No runtime database");
310  // Get Emc digitisation parameter container
311  fGeoPar = (PndEmcGeoPar*) db->getContainer("PndEmcGeoPar");
312  // Get Emc digitisation parameter container
313  fDigiPar = (PndEmcDigiPar*) db->getContainer("PndEmcDigiPar");
314  // Get Emc reconstruction parameter container
315  fRecoPar = (PndEmcRecoPar*) db->getContainer("PndEmcRecoPar");
316 }
Int_t run
Definition: autocutx.C:47
parameter set of Emc digitisation
Definition: PndEmcDigiPar.h:12
Parameter set for Emc Reco.
Definition: PndEmcRecoPar.h:12
void PndPersistencyTask::SetPersistency ( Bool_t  val = kTRUE)
inlineinherited

Definition at line 31 of file PndPersistencyTask.h.

References PndPersistencyTask::fPersistency, and val.

Referenced by barrelTrackFinder(), digi_complete(), digi_complete_newSTT(), digiOnly_complete(), PndBarrelTrackFinder::PndBarrelTrackFinder(), PndCATracking::PndCATracking(), PndDrcHitFinder::PndDrcHitFinder(), PndEmc2DLocMaxFinder::PndEmc2DLocMaxFinder(), PndEmcExpClusterSplitter::PndEmcExpClusterSplitter(), PndEmcFullDigiTask::PndEmcFullDigiTask(), PndEmcFWEndcapDigi::PndEmcFWEndcapDigi(), PndEmcFWEndcapTimebasedWaveforms::PndEmcFWEndcapTimebasedWaveforms(), PndEmcHitProducer::PndEmcHitProducer(), PndEmcHitsToWaveform::PndEmcHitsToWaveform(), PndEmcMakeBump::PndEmcMakeBump(), PndEmcMakeCluster::PndEmcMakeCluster(), PndEmcMakeClusterOnline::PndEmcMakeClusterOnline(), PndEmcMakeDigi::PndEmcMakeDigi(), PndEmcMakeRecoHit::PndEmcMakeRecoHit(), PndEmcMultiWaveformToCalibratedDigi::PndEmcMultiWaveformToCalibratedDigi(), PndEmcPhiBumpSplitter(), PndEmcTmpWaveformToDigi::PndEmcTmpWaveformToDigi(), PndEmcWaveformToCalibratedDigi::PndEmcWaveformToCalibratedDigi(), PndEmcWaveformToDigi::PndEmcWaveformToDigi(), PndFtofHitProducerIdeal::PndFtofHitProducerIdeal(), PndFtsCATracking::PndFtsCATracking(), PndFtsHitProducerIdeal::PndFtsHitProducerIdeal(), PndFtsHitProducerRealFast::PndFtsHitProducerRealFast(), PndFtsHitProducerRealFull::PndFtsHitProducerRealFull(), PndFtsHoughTrackerTask::PndFtsHoughTrackerTask(), PndGemDigitize::PndGemDigitize(), PndGemFindHits::PndGemFindHits(), PndIdealTrackFinder::PndIdealTrackFinder(), PndLmdPixelClusterTask::PndLmdPixelClusterTask(), PndLmdPixelHitProducerFast::PndLmdPixelHitProducerFast(), PndMdtClusterTask::PndMdtClusterTask(), PndMdtDigitization::PndMdtDigitization(), PndMdtHitProducerIdeal::PndMdtHitProducerIdeal(), PndMdtPointsToWaveform::PndMdtPointsToWaveform(), PndMdtTrkProducer::PndMdtTrkProducer(), PndMissingPzCleanerTask::PndMissingPzCleanerTask(), PndMvdGemTrackFinderOnHits::PndMvdGemTrackFinderOnHits(), PndMvdHitProducerIdeal::PndMvdHitProducerIdeal(), PndMvdPixelClusterTask::PndMvdPixelClusterTask(), PndMvdTimeWalkCorrTask::PndMvdTimeWalkCorrTask(), PndMvdToPix4ClusterTask::PndMvdToPix4ClusterTask(), PndRecoKalmanTask::PndRecoKalmanTask(), PndRecoKalmanTask2::PndRecoKalmanTask2(), PndRecoMultiKalmanTask::PndRecoMultiKalmanTask(), PndRecoMultiKalmanTask2::PndRecoMultiKalmanTask2(), PndRichHitFinder::PndRichHitFinder(), PndRichHitProducer::PndRichHitProducer(), PndSciTDigiTask::PndSciTDigiTask(), PndSciTHitProducerIdeal::PndSciTHitProducerIdeal(), PndSdsHitProducerIdeal::PndSdsHitProducerIdeal(), PndSdsHybridHitProducer::PndSdsHybridHitProducer(), PndSdsIdealClusterTask::PndSdsIdealClusterTask(), PndSdsIdealRecoTask::PndSdsIdealRecoTask(), PndSdsNoiseProducer::PndSdsNoiseProducer(), PndSdsPixelClusterTask::PndSdsPixelClusterTask(), PndSdsStripClusterTask::PndSdsStripClusterTask(), PndSdsStripHitProducer::PndSdsStripHitProducer(), PndSdsTimeWalkCorrTask::PndSdsTimeWalkCorrTask(), PndSttFindTracks::PndSttFindTracks(), PndSttHelixHitProducer::PndSttHelixHitProducer(), PndSttHitProducerIdeal::PndSttHitProducerIdeal(), PndSttHitProducerRealFast::PndSttHitProducerRealFast(), PndSttHitProducerRealFull::PndSttHitProducerRealFull(), PndSttMatchTracks::PndSttMatchTracks(), PndSttMvdGemTracking::PndSttMvdGemTracking(), PndTrackSmearTask::PndTrackSmearTask(), PndTrkTracking2::PndTrkTracking2(), reco(), reco_complete(), reco_complete_gf2(), reco_complete_newSTT(), reco_complete_sec(), recoideal_complete(), PndMvdClusterTask::SetPersistance(), PndMvdDigiTask::SetPersistance(), PndLmdPixelHitProducerFast::SetPersistance(), PndSdsHitProducerIdeal::SetPersistance(), PndSttMvdGemTracking::SetPersistenc(), PndMdtClusterTask::SetPersistence(), PndSttHelixHitProducer::SetPersistence(), PndMissingPzCleanerTask::SetPersistence(), PndFtsHitProducerRealFast::SetPersistence(), PndFtsHitProducerRealFull::SetPersistence(), PndSttHitProducerRealFull::SetPersistence(), PndSttHitProducerIdeal::SetPersistence(), PndSttHitProducerRealFast::SetPersistence(), PndFtsHitProducerIdeal::SetPersistence(), PndTrackSmearTask::SetPersistence(), PndSciTHitProducerIdeal::SetPersistence(), PndIdealTrackFinder::SetPersistence(), PndSttMatchTracks::SetPersistence(), PndSttFindTracks::SetPersistence(), PndFtsHoughTrackerTask::SetPersistence(), PndTrkTracking2::SetPersistence(), PndEmcMakeRecoHit::SetStorageOfData(), PndEmcFWEndcapDigi::SetStorageOfData(), PndEmcMakeClusterOnline::SetStorageOfData(), PndEmcFWEndcapTimebasedWaveforms::SetStorageOfData(), PndEmcMakeDigi::SetStorageOfData(), PndMdtPointsToWaveform::SetStorageOfData(), PndEmc2DLocMaxFinder::SetStorageOfData(), PndEmcMakeCluster::SetStorageOfData(), PndEmcHitsToWaveform::SetStorageOfData(), PndEmcMakeBump::SetStorageOfData(), PndEmcTmpWaveformToDigi::SetStorageOfData(), PndEmcWaveformToDigi::SetStorageOfData(), PndEmcWaveformToCalibratedDigi::SetStorageOfData(), PndEmcMultiWaveformToCalibratedDigi::SetStorageOfData(), PndEmcExpClusterSplitter::SetStorageOfData(), SetStorageOfData(), standard_tracking(), and PndEmcFullDigiTask::StoreDigi().

31 { fPersistency = val; }
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
void PndEmcPhiBumpSplitter::SetStorageOfData ( Bool_t  p = kTRUE)
inline

Definition at line 74 of file PndEmcPhiBumpSplitter.h.

References p, and PndPersistencyTask::SetPersistency().

74 {SetPersistency(p);};
Double_t p
Definition: anasim.C:58
void SetPersistency(Bool_t val=kTRUE)

Member Data Documentation

TClonesArray* PndEmcPhiBumpSplitter::fClusterArray
private

Definition at line 90 of file PndEmcPhiBumpSplitter.h.

Referenced by Exec(), and Init().

std::vector<Double_t> PndEmcPhiBumpSplitter::fClusterPosParam
private

Definition at line 101 of file PndEmcPhiBumpSplitter.h.

Referenced by Init(), and PndEmcPhiBumpSplitter().

TClonesArray* PndEmcPhiBumpSplitter::fDigiArray
private

Input array of PndEmcDigis

Definition at line 89 of file PndEmcPhiBumpSplitter.h.

Referenced by Exec(), and Init().

PndEmcDigiPar* PndEmcPhiBumpSplitter::fDigiPar
private

Definition at line 98 of file PndEmcPhiBumpSplitter.h.

Referenced by SetParContainers().

PndEmcGeoPar* PndEmcPhiBumpSplitter::fGeoPar
private

Definition at line 97 of file PndEmcPhiBumpSplitter.h.

Referenced by Init(), and SetParContainers().

TClonesArray* PndEmcPhiBumpSplitter::fPhiBumpArray
private

Output array of PndEmcBumps

Definition at line 93 of file PndEmcPhiBumpSplitter.h.

Referenced by AddPhiBump(), Exec(), and Init().

PndEmcRecoPar* PndEmcPhiBumpSplitter::fRecoPar
private

Definition at line 99 of file PndEmcPhiBumpSplitter.h.

Referenced by Init(), and SetParContainers().


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