Main Page | Modules | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members

CondensedNtpModule Class Reference

#include <CondensedNtpModule.h>

Inheritance diagram for CondensedNtpModule:

JobCModule List of all members.

Public Member Functions

 CondensedNtpModule ()
virtual ~CondensedNtpModule ()
void BeginJob ()
JobCResult Ana (const MomNavigator *mom)
const RegistryDefaultConfig () const
void Config (const Registry &r)
void Help ()
void EndJob ()

Private Member Functions

void FillUnRecoedMCInformation (TClonesArray *mcArray, TClonesArray *stdArray)
void ResetTreeVariables ()

Private Attributes

std::string fFileName
std::string fTreeName
TFile * fNtpFile
TTree * fNtuple
Int_t fDetector
Int_t fEndPlane
ANtpEventInfofEventInfo
ANtpHeaderInfofHeaderInfo
ANtpShowerInfofShowerInfo
ANtpTrackInfofTrackInfo
ANtpTruthInfofTruthInfo
ANtpInfoObjectFillerfInfoFiller
Int_t fStripInfoPresent
Int_t fDataType
Int_t fReconstructed

Constructor & Destructor Documentation

CondensedNtpModule::CondensedNtpModule  ) 
 

Definition at line 60 of file CondensedNtpModule.cxx.

References MSG.

00060                                        :
00061   fFileName("analysisNtuple.root"),
00062   fTreeName("analysisNtuple"),
00063   fNtpFile(0),
00064   fNtuple(0),
00065   fDetector(1),
00066   fEventInfo(0),
00067   fHeaderInfo(0),
00068   fShowerInfo(0),
00069   fTrackInfo(0),
00070   fTruthInfo(0),
00071   fStripInfoPresent(1),
00072   fDataType(1),
00073   fReconstructed(1)
00074 {
00075   MSG("JobC", Msg::kDebug) << "CondensedNtpModule::Constructor" << endl;
00076   
00077   fHeaderInfo = new ANtpHeaderInfo();
00078   fEventInfo = new ANtpEventInfo();
00079   fTrackInfo = new ANtpTrackInfo();
00080   fShowerInfo = new ANtpShowerInfo();
00081   fTruthInfo = new ANtpTruthInfo();
00082   
00083   //instantiate the ANtpInfoObjectFiller object
00084   fInfoFiller = new ANtpInfoObjectFiller();
00085         
00086 }

CondensedNtpModule::~CondensedNtpModule  )  [virtual]
 

Definition at line 89 of file CondensedNtpModule.cxx.

References MSG.

00090 {
00091   
00092   MSG("JobC", Msg::kDebug) << "CondensedNtpModule::Destructor" << endl;
00093   
00094   if(fHeaderInfo) delete fHeaderInfo;     
00095   if(fEventInfo) delete fEventInfo;       
00096   if(fShowerInfo) delete fShowerInfo;     
00097   if(fTrackInfo) delete fTrackInfo;       
00098   if(fTruthInfo) delete fTruthInfo;
00099     
00100 }


Member Function Documentation

JobCResult CondensedNtpModule::Ana const MomNavigator mom  )  [virtual]
 

Implement this for read only access to the MomNavigator

Reimplemented from JobCModule.

Definition at line 137 of file CondensedNtpModule.cxx.

References ANtpEventInfo::event, fDataType, fDetector, fEventInfo, ANtpInfoObjectFiller::FillEventInformation(), ANtpInfoObjectFiller::FillHeaderInformation(), ANtpInfoObjectFiller::FillMCTruthInformation(), ANtpInfoObjectFiller::FillShowerInformation(), ANtpInfoObjectFiller::FillTrackInformation(), FillUnRecoedMCInformation(), fInfoFiller, fNtuple, fReconstructed, fShowerInfo, fTrackInfo, fTruthInfo, ANtpRecoNtpManipulator::GetCosmicRayInfo(), ANtpEventManipulator::GetEvent(), ANtpRecoNtpManipulator::GetEventManipulator(), MomNavigator::GetFragment(), ANtpRecoNtpManipulator::GetMCArray(), ANtpRecoNtpManipulator::GetMCManipulator(), ANtpMCManipulator::GetNtpMCStdHep(), ANtpMCManipulator::GetNtpMCTruth(), ANtpMCManipulator::GetNtpTHEvent(), ANtpEventManipulator::GetPrimaryShower(), ANtpEventManipulator::GetPrimaryTrack(), ANtpRecoNtpManipulator::GetRun(), ANtpRecoNtpManipulator::GetSnarl(), ANtpRecoNtpManipulator::GetSnarlEventSummary(), ANtpRecoNtpManipulator::GetStdHepArray(), ANtpEventManipulator::GetStrip(), ANtpRecoNtpManipulator::GetSubRun(), MSG, NtpTHEvent::neumc, NtpTHEvent::neustdhep, NtpSREventSummary::nevent, ResetTreeVariables(), ANtpEventManipulator::SetEventInSnarl(), JobCResult::SetFailed(), ANtpRecoNtpManipulator::SetPrimaryShowerCriteria(), ANtpRecoNtpManipulator::SetPrimaryTrackCriteria(), and JobCResult::SetWarning().

00138 {
00139   MSG("CondensedNtpModule", Msg::kDebug) << "start ana" << endl;
00140   
00141   JobCResult result(JobCResult::kPassed);
00142   
00143   //get the records from MOM
00144   assert(mom);
00145   
00146   NtpSRRecord *record = dynamic_cast<NtpSRRecord *>(mom->GetFragment("NtpSRRecord"));
00147   NtpStRecord *stRecord = dynamic_cast<NtpStRecord *>(mom->GetFragment("NtpStRecord"));
00148   
00149   if(!record && !stRecord){
00150     MSG("CondensedNtpModuleAtm", Msg::kWarning) << "Could not get NtpSR or NtpSt Record from MOM" << endl;
00151     result.SetWarning().SetFailed();
00152     return result;
00153   }//end if no NtpSRRecord
00154   
00155   //if there is no NtpMCRecord, then mom returns a null pointer
00156   NtpMCRecord *mcRecord = 0;
00157   mcRecord = dynamic_cast<NtpMCRecord *>(mom->GetFragment("NtpMCRecord"));
00158   if(!mcRecord&&fDataType&&!stRecord){
00159     MSG("CondensedNtpModule", Msg::kWarning) << "Could not get NtpMCRecord from MOM " 
00160                                              << " - hopefully you werent expecting one " << endl;
00161     result.SetWarning().SetFailed();
00162     return result;
00163   }
00164   
00165   NtpTHRecord *thRecord = 0;
00166   thRecord = dynamic_cast<NtpTHRecord *>(mom->GetFragment("NtpTHRecord"));
00167   
00168   //fill the header information
00169   DetectorType::Detector_t det = DetectorType::kFar;
00170   if(fDetector<1) det = DetectorType::kNear;
00171   //get a sim flag
00172   SimFlag::SimFlag_t dataType = SimFlag::kData;
00173   if(fDataType) dataType = SimFlag::kMC;
00174   
00175   //instantiate a ANtpRecoNtpManipulator object to help you get the info you want
00176   ANtpRecoNtpManipulator *ntpManipulator = 0;
00177   if(record) ntpManipulator = new ANtpRecoNtpManipulator(record, mcRecord, thRecord);
00178   else if(stRecord) ntpManipulator = new ANtpRecoNtpManipulator(stRecord);
00179 
00180   fInfoFiller->FillHeaderInformation(ntpManipulator->GetRun(), ntpManipulator->GetSubRun(),
00181                                      ntpManipulator->GetSnarl(), ntpManipulator->GetSnarlEventSummary(), 
00182                                      ntpManipulator->GetCosmicRayInfo(), det, dataType, fHeaderInfo);  
00183 
00184   //set up which flags you want to use to determine the primary shower or track
00185   //a value of 0 for a flag means it will not be used
00186   ntpManipulator->SetPrimaryTrackCriteria(0,1,0); // nplanes, length, total pulse height
00187   ntpManipulator->SetPrimaryShowerCriteria(0,1); // nplanes, total pulse height
00188   
00189   //mc data because fDataType = 1, but no reconstructed event because 
00190   //record doesnt exist.  still want to have the mc info in the tree though 
00191   //to understand events that werent reconstructed
00192   if(mcRecord&&!record&&fDataType&&!stRecord){
00193     
00194     MSG("CondensedNtpModule", Msg::kDebug) << "why am i here?" << endl;
00195     
00196     fReconstructed = 0;
00197     
00198     //flag the run, subRun values to indicated that 
00199     //this snarl was not reconstructed.  keep track of 
00200     //how many times that happens
00201     FillUnRecoedMCInformation(ntpManipulator->GetMCArray(), ntpManipulator->GetStdHepArray());
00202     fNtuple->Fill();
00203     ResetTreeVariables();
00204     
00205     return result;              
00206   }
00207   
00208   fReconstructed = 1;
00209   
00210   //-----------------------------------------------------------------------------------
00211   //here is where you start looping over the events in the ntuple to fill your analysis tree
00212   
00213   //set up some pointers to typical NtpSR objects
00214   NtpSRTrack *ntpTrack = 0;
00215   NtpSRStrip *ntpStrip = 0;
00216   NtpSREvent *ntpEvent = 0;
00217   NtpSRShower *ntpShower = 0;
00218   NtpMCTruth *ntpMCTruth = 0;
00219   NtpMCStdHep *ntpMCStdHep = 0;
00220   NtpTHEvent *ntpTHEvent = 0;
00221   //NtpTHStrip *ntpTHStrip = 0;
00222   //NtpTHTrack *ntpTHTrack = 0;
00223   //NtpTHShower *ntpTHShower = 0;
00224   
00225   //loop over the events in the snarl - it is also possible to just loop over 
00226   //showers, tracks or slices in the record using a loop like
00227   //for(Int_t i = 0; i < record->evthdr.ntrack; ++i){
00228   //  ntpTrack = ntpManipulator->GetSnarlTrack(i);
00229   //}
00230   //
00231   //you can get the primary track/shower from the snarl using 
00232   //ntpManipulator->GetSnarlPrimaryTrack()/GetSnarlPrimaryShower()
00233   
00234   for(Int_t i = 0; i < ntpManipulator->GetSnarlEventSummary().nevent; ++i){
00235     
00236     MSG("CondensedNtpModule", Msg::kDebug) << "on event " << i << endl;
00237     
00238     fEventInfo->event = i;
00239     
00240     //get event i.  the call to NtpHelper::GetEventInSnarl sets the 
00241     //NtpSREvent data member in that object to the current event so that
00242     //you can later call for the primary track, shower, etc.
00243     ntpManipulator->GetEventManipulator()->SetEventInSnarl(i);
00244     ntpEvent = ntpManipulator->GetEventManipulator()->GetEvent();               
00245     if(ntpEvent) fInfoFiller->FillEventInformation(ntpEvent, fEventInfo);
00246     
00247     //get the primary track for the event - if no track is present it 
00248     //returns 0
00249     ntpTrack = ntpManipulator->GetEventManipulator()->GetPrimaryTrack();                
00250     if(ntpTrack) fInfoFiller->FillTrackInformation(ntpTrack, fTrackInfo);
00251     
00252     //get the primary shower for the event - if no shower is present it 
00253     //returns 0
00254     ntpShower = ntpManipulator->GetEventManipulator()->GetPrimaryShower();              
00255     if(ntpShower) fInfoFiller->FillShowerInformation(ntpShower, fShowerInfo);
00256     
00257     //why not get a strip as well - if there is a strip in the event, the 
00258     //following line will get it
00259     if(fStripInfoPresent) ntpStrip = ntpManipulator->GetEventManipulator()->GetStrip(0);
00260     
00261     MSG("CondensedNtpModule", Msg::kDebug) << "CondensedNtpModule::Ana mcRecord="<<mcRecord << endl;
00262     MSG("CondensedNtpModule", Msg::kDebug) << "CondensedNtpModule::Ana thRecord="<<thRecord << endl;
00263     //check to see if there is MC information around and get the NtpMCTruth and
00264     //NtpTHEvent objects
00265 
00266     //the GetNtpTHEvent method takes as its argument the same index as the 
00267     //GetEventInSnarl method, ie the index of the current event.  you can 
00268     //then use the NtpTHEvent object returned to figure out which index to 
00269     //give the GetNtpMCTruth and GetNtpMCStdHep methods to return the appropriate 
00270     //truth information corresponding to this event.  what is going on behind the 
00271     //scenes is that there are arrays of NtpMCTruth and NtpMCStdHep objects in 
00272     //the NtpMCRecord but to figure out which entry is the appropriate one you need 
00273     //to have the NtpTHRecord array of thevt which is filled in the same order as 
00274     //the evt array in the NtpSRRecord.
00275     ntpTHEvent = ntpManipulator->GetMCManipulator()->GetNtpTHEvent(i);
00276     if(ntpTHEvent){
00277       ntpMCTruth = ntpManipulator->GetMCManipulator()->GetNtpMCTruth(ntpTHEvent->neumc);        
00278       ntpMCStdHep = ntpManipulator->GetMCManipulator()->GetNtpMCStdHep(ntpTHEvent->neustdhep);  
00279     }
00280     
00281     //can also do the above using NtpTHStrips, NtpTHShowers or NtpTHTracks.  get the index 
00282     //of the strip, track or shower in the NtpSRRecord stp, trk or shw arrays first, then 
00283     //get the NtpTHTrack/Shower object using the NtpHelper.  have to uncomment the variable 
00284     //declarations in lines 223-225 to use this
00285     //uncomment out the disired lines below
00286     /*
00287     
00288       if(ntpStrip) ntpTHStrip = ntpManipulator->GetMCManipulator()->GetNtpTHStrip(ntpStrip->index);
00289       if(ntpTHStrip){
00290       ntpMCTruth = ntpManipulator->GetMCManipulator()->GetNtpMCTruth(ntpTHStrip->neumc);        
00291       ntpMCStdHep = ntpManipulator->GetMCManipulator()->GetNtpMCStdHep(ntpTHStrip->neustdhep);  
00292       }
00293     
00294       if(ntpTrack) ntpTHTrack = ntpManipulator->GetMCManipulator()->GetNtpTHTrack((Int_t)ntpTrack->index);
00295       if(ntpTHTrack){
00296       ntpMCTruth = ntpManipulator->GetMCManipulator()->GetNtpMCTruth(ntpTHTrack->neumc);        
00297       ntpMCStdHep = ntpManipulator->GetMCManipulator()->GetNtpMCStdHep(ntpTHTrack->neustdhep);  
00298       }
00299       if(ntpShower) ntpTHShower = ntpManipulator->GetMCManipulator()->GetNtpTHShower((Int_t)ntpShower->index);   
00300       if(ntpTHShower){
00301       ntpMCTruth = ntpManipulator->GetMCManipulator()->GetNtpMCTruth(ntpTHShower->neumc);       
00302       ntpMCStdHep = ntpManipulator->GetMCManipulator()->GetNtpMCStdHep(ntpTHShower->neustdhep); 
00303       }
00304     */
00305       
00306     //got the NtpMCTruth object now, so fill the desired info
00307     MSG("CondensedNtpModule", Msg::kDebug) << "CondensedNtpModule::Ana ntpMCTruth="<<ntpMCTruth<< endl;
00308     if(ntpMCTruth){
00309       fInfoFiller->FillMCTruthInformation(ntpMCTruth, fTruthInfo);
00310     }                           
00311   
00312     
00313     fNtuple->Fill();
00314     ResetTreeVariables();
00315   }//end loop over events for this snarl
00316   //-----------------------------------------------------------------------------------
00317 
00318   delete ntpManipulator;
00319 
00320   return result;
00321 }

void CondensedNtpModule::BeginJob  )  [virtual]
 

Implement for notification of begin of job

Reimplemented from JobCModule.

Definition at line 103 of file CondensedNtpModule.cxx.

References fDataType, fEventInfo, fFileName, fHeaderInfo, fNtpFile, fNtuple, fShowerInfo, fTrackInfo, fTreeName, fTruthInfo, and MSG.

00104 {
00105   // save the current working directory
00106   TDirectory* savedir  = gDirectory;  
00107   
00108   MSG("CondensedNtpModule", Msg::kDebug) << fFileName.c_str() << " " << fTreeName.c_str() << endl;
00109   
00110   //create the new TFile for holding the electronics tree
00111   // this changes the value of gDirectoy
00112   fNtpFile = new TFile(fFileName.c_str(),"RECREATE");  
00113   
00114   // create TTree, these will attach themselves to the current 
00115   //working directory   
00116   fNtuple = new TTree(fTreeName.c_str(), "Analysis Tree");
00117   
00118   //the branches are all objects of type ANtp*Info*
00119   fNtuple->Branch("header.", "ANtpHeaderInfo", &fHeaderInfo, 64000, 2);
00120   fNtuple->Branch("event.", "ANtpEventInfo", &fEventInfo, 64000, 2);
00121   fNtuple->Branch("shower.", "ANtpShowerInfo", &fShowerInfo, 64000, 2);
00122   fNtuple->Branch("track.", "ANtpTrackInfo", &fTrackInfo, 64000, 2);
00123   if(fDataType>0){
00124     fNtuple->Branch("truth.", "ANtpTruthInfo", &fTruthInfo, 64000, 2);  
00125   } 
00126   //---------------------------------------------------------------------------------  
00127   
00128   // change back to current working directory before leaving constructor
00129   savedir->cd();   
00130   
00131   return;
00132 }

void CondensedNtpModule::Config const Registry r  )  [virtual]
 

Return the actual configuration. If your module directly pulls its configuration from the fConfig Registry, you don't need to override this. Override if you have local config variables.

Reimplemented from JobCModule.

Definition at line 383 of file CondensedNtpModule.cxx.

References fDataType, fDetector, fEndPlane, fFileName, fStripInfoPresent, fTreeName, and Registry::Get().

00384 {
00385   int         tmpb;  // a temp bool. See comment under DefaultConfig...
00386   int         tmpi;  // a temp int.
00387   double      tmpd;  // a temp double.
00388   const char* tmps;  // a temp string
00389   
00390   if (r.Get("FileName",           tmps)) fFileName              =  tmps;
00391   if (r.Get("TreeName",           tmps)) fTreeName              =  tmps;
00392   if (r.Get("EndPlane",           tmpi)) fEndPlane              =  tmpi;
00393   if (r.Get("Detector",           tmpi)) fDetector              =  tmpi;
00394   if (r.Get("StripInfoPresent",   tmpi)) fStripInfoPresent      =  tmpi;
00395   if (r.Get("DataType",           tmpi)) fDataType              =  tmpi;
00396   
00397   //quiet compiler warnings
00398   tmpb = 0;
00399   tmpd = 0.;
00400   
00401   return;
00402 }

const Registry & CondensedNtpModule::DefaultConfig  )  const [virtual]
 

Get the default configuration registry. This should normally be overridden. One useful idiom is to implement it like:

const Registry& MyModule::DefaultConfig() const { static Registry cfg; // never is destroyed if (cfg.Size()) return cfg; // already filled it // set defaults: cfg.Set("TheAnswer",42); cfg.Set("Units","unknown"); return cfg; }

Reimplemented from JobCModule.

Definition at line 357 of file CondensedNtpModule.cxx.

References Registry::LockValues(), Registry::Set(), and Registry::UnLockValues().

00358 {
00359   
00360   int itrue = 1;  // work around for lack of bool in registry
00361   int ifalse = 0; // work around for lack of bool in registry
00362   
00363   static Registry r;
00364   
00365   r.UnLockValues();
00366   
00367   r.Set("FileName",           "analysisTree.root");
00368   r.Set("TreeName",           "analysisTree");
00369   r.Set("EndPlane",           485);
00370   r.Set("Detector",           1);
00371   r.Set("StripInfoPresent",   1);
00372   r.Set("DataType",           1);
00373   
00374   r.LockValues();
00375   
00376   //quiet compiler warnings
00377   itrue = ifalse;
00378   
00379   return r;
00380 }

void CondensedNtpModule::EndJob  )  [virtual]
 

Implement for notification of end of job

Reimplemented from JobCModule.

Definition at line 334 of file CondensedNtpModule.cxx.

References fNtpFile, fNtuple, and MSG.

00335 {
00336   
00337   MSG("CondensedNtpModule", Msg::kInfo) << "start end job method" << endl;      
00338   
00339   //save current working directory
00340   TDirectory *savedir = gDirectory;
00341   
00342   //cd to the "directory" of the file
00343   fNtpFile->cd();
00344 
00345   fNtuple->Write();
00346   
00347   //go back to original working directory
00348   savedir->cd();  
00349   
00350   fNtpFile->Write();
00351   fNtpFile->Close();
00352   
00353   return;
00354 }

void CondensedNtpModule::FillUnRecoedMCInformation TClonesArray *  mcArray,
TClonesArray *  stdArray
[private]
 

Definition at line 405 of file CondensedNtpModule.cxx.

References ANtpEventInfo::event, fEventInfo, ANtpInfoObjectFiller::FillMCTruthInformation(), fInfoFiller, fTruthInfo, and NtpMCTruth::stdhep.

Referenced by Ana().

00406 {
00407   NtpMCTruth *ntpMCTruth = 0;
00408   NtpMCStdHep *ntpMCStdHep = 0;
00409   
00410   //get the number of entries in the arrays
00411   Int_t numMCEntries = mcArray->GetLast()+1;
00412   Int_t stdHepLow = 0;
00413   Int_t stdHepHigh = 0;
00414   
00415   //loop over the individual events - only one in the far detector, 
00416   //many in the near
00417   for(Int_t i = 0; i < numMCEntries; ++i){
00418     
00419     fEventInfo->event = i;
00420     
00421     ntpMCTruth = dynamic_cast<NtpMCTruth *>(mcArray->At(i));
00422     if(ntpMCTruth) fInfoFiller->FillMCTruthInformation(ntpMCTruth, fTruthInfo);
00423     
00424     stdHepLow = ntpMCTruth->stdhep[0];
00425     stdHepHigh = ntpMCTruth->stdhep[1];
00426     
00427     //loop over the stdhep entries for the current mc event
00428     for(Int_t j = stdHepLow; j < stdHepHigh+1; ++j){
00429       ntpMCStdHep = dynamic_cast<NtpMCStdHep *>(stdArray->At(j));
00430       
00431       //add the code for what you want to fill here 
00432     }                   
00433     
00434   }//end loop over mc events
00435   
00436   
00437   return;
00438 }

void CondensedNtpModule::Help  )  [virtual]
 

Implement to spew some useful help to cout

Reimplemented from JobCModule.

Definition at line 324 of file CondensedNtpModule.cxx.

References MSG.

00325 {
00326   MSG("JobC", Msg::kInfo) 
00327     << "NearElectronicsCheck::Help\n"
00328     <<"CondensedNtpModule is a module which demultiplexes events "
00329     <<"in the far detector."
00330     << endl;
00331 }

void CondensedNtpModule::ResetTreeVariables  )  [private]
 

Definition at line 441 of file CondensedNtpModule.cxx.

References fEventInfo, fShowerInfo, fTrackInfo, fTruthInfo, ANtpTruthInfo::Reset(), ANtpTrackInfo::Reset(), ANtpShowerInfo::Reset(), and ANtpEventInfo::Reset().

Referenced by Ana().

00442 {
00443     
00444   fEventInfo->Reset();
00445   fShowerInfo->Reset();
00446   fTrackInfo->Reset();
00447   fTruthInfo->Reset();
00448   
00449   return;
00450 }


Member Data Documentation

Int_t CondensedNtpModule::fDataType [private]
 

Definition at line 66 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), and Config().

Int_t CondensedNtpModule::fDetector [private]
 

Definition at line 54 of file CondensedNtpModule.h.

Referenced by Ana(), and Config().

Int_t CondensedNtpModule::fEndPlane [private]
 

Definition at line 55 of file CondensedNtpModule.h.

Referenced by Config().

ANtpEventInfo* CondensedNtpModule::fEventInfo [private]
 

Definition at line 57 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), FillUnRecoedMCInformation(), and ResetTreeVariables().

std::string CondensedNtpModule::fFileName [private]
 

Definition at line 50 of file CondensedNtpModule.h.

Referenced by BeginJob(), and Config().

ANtpHeaderInfo* CondensedNtpModule::fHeaderInfo [private]
 

Definition at line 58 of file CondensedNtpModule.h.

Referenced by BeginJob().

ANtpInfoObjectFiller* CondensedNtpModule::fInfoFiller [private]
 

Definition at line 63 of file CondensedNtpModule.h.

Referenced by Ana(), and FillUnRecoedMCInformation().

TFile* CondensedNtpModule::fNtpFile [private]
 

Definition at line 52 of file CondensedNtpModule.h.

Referenced by BeginJob(), and EndJob().

TTree* CondensedNtpModule::fNtuple [private]
 

Definition at line 53 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), and EndJob().

Int_t CondensedNtpModule::fReconstructed [private]
 

Definition at line 67 of file CondensedNtpModule.h.

Referenced by Ana().

ANtpShowerInfo* CondensedNtpModule::fShowerInfo [private]
 

Definition at line 59 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), and ResetTreeVariables().

Int_t CondensedNtpModule::fStripInfoPresent [private]
 

Definition at line 65 of file CondensedNtpModule.h.

Referenced by Config().

ANtpTrackInfo* CondensedNtpModule::fTrackInfo [private]
 

Definition at line 60 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), and ResetTreeVariables().

std::string CondensedNtpModule::fTreeName [private]
 

Definition at line 51 of file CondensedNtpModule.h.

Referenced by BeginJob(), and Config().

ANtpTruthInfo* CondensedNtpModule::fTruthInfo [private]
 

Definition at line 61 of file CondensedNtpModule.h.

Referenced by Ana(), BeginJob(), FillUnRecoedMCInformation(), and ResetTreeVariables().


The documentation for this class was generated from the following files:
Generated on Thu Nov 1 15:56:15 2007 for loon by  doxygen 1.3.9.1