#include <FarDetShieldPlankListModule.h>
Inheritance diagram for FarDetShieldPlankListModule:

Public Member Functions | |
| FarDetShieldPlankListModule () | |
| ~FarDetShieldPlankListModule () | |
| void | BeginJob () |
| JobCResult | Reco (MomNavigator *mom) |
| JobCResult | Ana (const MomNavigator *mom) |
| const Registry & | DefaultConfig () const |
| void | Config (const Registry &r) |
| void | EndJob () |
Private Attributes | |
| TString | fListIn |
| TString | fListOut |
| Int_t | run |
| Int_t | snarl |
| Int_t | timeframe |
| Int_t | time |
| Int_t | date |
| Int_t | section |
| Int_t | subsection |
| Int_t | plane |
| Int_t | Nstrips |
| Int_t | Ndigits |
| Int_t | Nerrors |
| Double_t | Xpos |
| Double_t | Ypos |
| Int_t | crateN |
| Int_t | varcN |
| Int_t | vmmN |
| Int_t | vaadcN |
| Int_t | vachipN |
| Int_t | vachannelN |
| Double_t | TrawN |
| Double_t | TcalN |
| Double_t | QrawN |
| Double_t | QcalN |
| Double_t | ZposN |
| Int_t | crateS |
| Int_t | varcS |
| Int_t | vmmS |
| Int_t | vaadcS |
| Int_t | vachipS |
| Int_t | vachannelS |
| Double_t | TrawS |
| Double_t | TcalS |
| Double_t | QrawS |
| Double_t | QcalS |
| Double_t | ZposS |
| TFile * | fShieldFile |
| TTree * | fShieldTree |
|
|
Definition at line 31 of file FarDetShieldPlankListModule.cxx. 00031 : 00032 fListIn("canddigitlist"), 00033 fListOut("FarDetShieldPlankListHandle"), 00034 fShieldFile(0), 00035 fShieldTree(0) 00036 { 00037 00038 }
|
|
|
Definition at line 40 of file FarDetShieldPlankListModule.cxx. 00041 {
00042
00043 }
|
|
|
Implement this for read only access to the MomNavigator Reimplemented from JobCModule. Definition at line 91 of file FarDetShieldPlankListModule.cxx. References crateN, crateS, date, digit(), CandRecord::FindCandHandle(), RawRecord::FindRawBlock(), fShieldFile, fShieldTree, CandRecord::GetCandHeader(), CandDigitHandle::GetChannelId(), FarDetShieldPlankHandle::GetCharge(), RawChannelId::GetCrate(), CandHandle::GetDaughterIterator(), PlexSEIdAltL::GetEnd(), MomNavigator::GetFragment(), FarDetShieldPlankHandle::GetGeomErrors(), CandHandle::GetNDaughters(), FarDetShieldPlankHandle::GetNStrips(), FarDetShieldPlankHandle::GetPlane(), CandDigitHandle::GetPlexSEIdAltL(), RawRecord::GetRawHeader(), CandHeader::GetRun(), RawDaqHeader::GetRun(), FarDetShieldPlankHandle::GetSection(), CandHeader::GetSnarl(), RawDaqSnarlHeader::GetSnarl(), FarDetShieldPlankHandle::GetSubSection(), FarDetShieldPlankHandle::GetTime(), RawDaqHeader::GetTimeFrameNum(), RawSnarlHeaderBlock::GetTriggerTime(), RawChannelId::GetVaAdcSel(), RawChannelId::GetVaChannel(), RawChannelId::GetVaChip(), RawChannelId::GetVarcId(), RawChannelId::GetVmm(), FarDetShieldPlankHandle::GetX(), FarDetShieldPlankHandle::GetY(), FarDetShieldPlankHandle::GetZ(), MSG, Ndigits, Nerrors, Nstrips, plane, QcalN, QcalS, QrawN, QrawS, run, section, JobCResult::SetFailed(), snarl, subsection, TcalN, TcalS, time, timeframe, TrawN, TrawS, vaadcN, vaadcS, vachannelN, vachannelS, vachipN, vachipS, varcN, varcS, vmmN, vmmS, Xpos, Ypos, ZposN, and ZposS. 00092 {
00093 MSG("FarDetShieldPlank",Msg::kDebug) << " *** FarDetShieldPlankListModule::Ana(...) *** " << endl;
00094
00095 JobCResult result(JobCResult::kPassed);
00096
00097 run=-1; snarl=-1; timeframe=-1; time=-1; date=-1;
00098
00099 const RawRecord *rawrec = dynamic_cast<const RawRecord*>(mom->GetFragment("RawRecord"));
00100 if(rawrec == 0){
00101 MSG("FarDetShieldPlank", Msg::kWarning) << " *** FAILED TO FIND RAWRECORD *** " << endl;
00102 return result.SetFailed();
00103 }
00104
00105 const RawDaqSnarlHeader* rawheader = dynamic_cast<const RawDaqSnarlHeader*>(rawrec->GetRawHeader());
00106 if(rawheader){
00107 if(run<0) run = rawheader->GetRun();
00108 if(snarl<0) snarl = rawheader->GetSnarl();
00109 if(timeframe<0) timeframe = rawheader->GetTimeFrameNum();
00110 }
00111
00112 const RawSnarlHeaderBlock* rdb = dynamic_cast<const RawSnarlHeaderBlock*>(rawrec->FindRawBlock("RawSnarlHeaderBlock"));
00113 if(rdb){
00114 date = (((VldTimeStamp)(rdb->GetTriggerTime())).GetSec()-1059696000)/(3600*24);
00115 time = (((VldTimeStamp)(rdb->GetTriggerTime())).GetSec()-1059696000)%(3600*24);
00116 }
00117
00118 const CandRecord *const candrec = dynamic_cast<const CandRecord *>(mom->GetFragment("CandRecord", "PrimaryCandidateRecord"));
00119 if(candrec == 0){
00120 MSG("FarDetShieldPlank", Msg::kWarning) << " *** FAILED TO FIND CANDRECORD *** " << endl;
00121 return result.SetFailed();
00122 }
00123
00124 const CandHeader* candheader = dynamic_cast<const CandHeader*>(candrec->GetCandHeader());
00125 if(candheader){
00126 if(run<0) run = candheader->GetRun();
00127 if(snarl<0) snarl = candheader->GetSnarl();
00128 }
00129
00130 const FarDetShieldPlankListHandle* shldlist = dynamic_cast<const FarDetShieldPlankListHandle*>(candrec->FindCandHandle("FarDetShieldPlankListHandle"));
00131 if(shldlist==0){
00132 return result.SetFailed();
00133 }
00134
00135 if(!fShieldFile){
00136 // TString mystring("results/shieldplank");
00137 // mystring.Append("."); mystring+=run;
00138 // mystring.Append(".root");
00139 TString mystring("shieldplank.root");
00140 TDirectory* tmpd = gDirectory;
00141 fShieldFile = new TFile(mystring.Data(),"RECREATE");
00142 fShieldTree = new TTree("ShieldTree","ShieldTree");
00143 fShieldTree->SetAutoSave(100);
00144
00145 fShieldTree->Branch("run",&run,"run/I");
00146 fShieldTree->Branch("snarl",&snarl,"snarl/I");
00147 fShieldTree->Branch("timeframe",&timeframe,"timeframe/I");
00148 fShieldTree->Branch("time",&time,"time/I");
00149 fShieldTree->Branch("date",&date,"date/I");
00150 fShieldTree->Branch("section",§ion,"section/I");
00151 fShieldTree->Branch("subsection",&subsection,"subsection/I");
00152 fShieldTree->Branch("plane",&plane,"plane/I");
00153 fShieldTree->Branch("Nerrors",&Nerrors,"Nerrors/I");
00154 fShieldTree->Branch("Nstrips",&Nstrips,"Nstrips/I");
00155 fShieldTree->Branch("Ndigits",&Ndigits,"Ndigits/I");
00156 fShieldTree->Branch("Xpos",&Xpos,"Xpos/D");
00157 fShieldTree->Branch("Ypos",&Ypos,"Ypos/D");
00158
00159 fShieldTree->Branch("crateN",&crateN,"crateN/I");
00160 fShieldTree->Branch("varcN",&varcN,"varcN/I");
00161 fShieldTree->Branch("vmmN",&vmmN,"vmmN/I");
00162 fShieldTree->Branch("vaadcN",&vaadcN,"vaadcN/I");
00163 fShieldTree->Branch("vachipN",&vachipN,"vachipN/I");
00164 fShieldTree->Branch("vachannelN",&vachannelN,"vachannelN/I");
00165 fShieldTree->Branch("TrawN",&TrawN,"TrawN/D");
00166 fShieldTree->Branch("TcalN",&TcalN,"TcalN/D");
00167 fShieldTree->Branch("QrawN",&QrawN,"QrawN/D");
00168 fShieldTree->Branch("QcalN",&QcalN,"QcalN/D");
00169 fShieldTree->Branch("ZposN",&ZposN,"ZposN/D");
00170
00171 fShieldTree->Branch("crateS",&crateS,"crateS/I");
00172 fShieldTree->Branch("varcS",&varcS,"varcS/I");
00173 fShieldTree->Branch("vmmS",&vmmS,"vmmS/I");
00174 fShieldTree->Branch("vaadcS",&vaadcS,"vaadcS/I");
00175 fShieldTree->Branch("vachipS",&vachipS,"vachipS/I");
00176 fShieldTree->Branch("vachannelS",&vachannelS,"vachannelS/I");
00177 fShieldTree->Branch("TrawS",&TrawS,"TrawS/D");
00178 fShieldTree->Branch("TcalS",&TcalS,"TcalS/D");
00179 fShieldTree->Branch("QrawS",&QrawS,"QrawS/D");
00180 fShieldTree->Branch("QcalS",&QcalS,"QcalS/D");
00181 fShieldTree->Branch("ZposS",&ZposS,"ZposS/D");
00182
00183 gDirectory = tmpd;
00184 }
00185
00186 TIter shlditr(shldlist->GetDaughterIterator());
00187 while(FarDetShieldPlankHandle* shld = dynamic_cast<FarDetShieldPlankHandle*>(shlditr())){
00188
00189 section = shld->GetSection();
00190 subsection = shld->GetSubSection();
00191 plane = shld->GetPlane();
00192 Nstrips = shld->GetNStrips();
00193 Ndigits = shld->GetNDaughters();
00194 Nerrors = shld->GetGeomErrors();
00195 Xpos = shld->GetX();
00196 Ypos = shld->GetY();
00197
00198 crateN = -1;
00199 varcN = -1;
00200 vmmN = -1;
00201 vaadcN = -1;
00202 vachipN = -1;
00203 vachannelN = -1;
00204 TrawN = 0.0;
00205 TcalS = 0.0;
00206 QrawN = 0.0;
00207 QcalN = 0.0;
00208 ZposN = 0.0;
00209
00210 crateS = -1;
00211 varcS = -1;
00212 vmmS = -1;
00213 vaadcS = -1;
00214 vachipS = -1;
00215 vachannelS = -1;
00216 TrawS = 0.0;
00217 TcalS = 0.0;
00218 QrawS = 0.0;
00219 QcalS = 0.0;
00220 ZposS = 0.0;
00221
00222 TIter digitr(shld->GetDaughterIterator());
00223 while(CandDigitHandle* digit = dynamic_cast<CandDigitHandle*>(digitr())){
00224
00225 if(digit->GetPlexSEIdAltL().GetEnd()==StripEnd::kNegative){
00226 TrawN = shld->GetTime(StripEnd::kNegative,CalTimeType::kNone);
00227 TcalN = shld->GetTime(StripEnd::kNegative,CalTimeType::kT0);
00228 QrawN = shld->GetCharge(StripEnd::kNegative,CalDigitType::kNone);
00229 QcalN = shld->GetCharge(StripEnd::kNegative,CalDigitType::kPE);
00230 ZposN = shld->GetZ(StripEnd::kNegative);
00231
00232 RawChannelId rawch = digit->GetChannelId();
00233 crateN=rawch.GetCrate();
00234 varcN=rawch.GetVarcId();
00235 vmmN=rawch.GetVmm();
00236 vaadcN=rawch.GetVaAdcSel();
00237 vachipN=rawch.GetVaChip();
00238 vachannelN=rawch.GetVaChannel();
00239 }
00240
00241 if(digit->GetPlexSEIdAltL().GetEnd()==StripEnd::kPositive){
00242 TrawS = shld->GetTime(StripEnd::kPositive,CalTimeType::kNone);
00243 TcalS = shld->GetTime(StripEnd::kPositive,CalTimeType::kT0);
00244 QrawS = shld->GetCharge(StripEnd::kPositive,CalDigitType::kNone);
00245 QcalS = shld->GetCharge(StripEnd::kPositive,CalDigitType::kPE);
00246 ZposS = shld->GetZ(StripEnd::kPositive);
00247
00248 RawChannelId rawch = digit->GetChannelId();
00249 crateS=rawch.GetCrate();
00250 varcS=rawch.GetVarcId();
00251 vmmS=rawch.GetVmm();
00252 vaadcS=rawch.GetVaAdcSel();
00253 vachipS=rawch.GetVaChip();
00254 vachannelS=rawch.GetVaChannel();
00255 }
00256
00257 }
00258
00259 TDirectory* tmpd = gDirectory;
00260 fShieldFile->cd();
00261 fShieldTree->Fill();
00262 gDirectory = tmpd;
00263 }
00264
00265 return result;
00266 }
|
|
|
Implement for notification of begin of job Reimplemented from JobCModule. Definition at line 45 of file FarDetShieldPlankListModule.cxx. References MSG. 00046 {
00047 MSG("FarDetShieldPlank",Msg::kDebug) << " *** FarDetShieldPlankListModule::BeginJob() *** " << endl;
00048
00049 /*
00050 AlgFactory &af = AlgFactory::GetInstance();
00051 af.Register("AlgFarDetShieldPlankList","default","libFarDetShieldPlank.so","AlgConfig");
00052 af.Register("AlgFarDetShieldPlank","default","libFarDetShieldPlank.so","AlgConfig");
00053 */
00054 }
|
|
|
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 296 of file FarDetShieldPlankListModule.cxx. References fListIn, fListOut, Registry::Get(), AlgHandle::GetAlgConfig(), AlgFactory::GetAlgHandle(), AlgFactory::GetInstance(), Registry::LockValues(), MSG, and Registry::UnLockValues(). 00297 {
00298 MSG("FarDetShieldPlank",Msg::kDebug) << " *** FarDetShieldPlankListModule::Config( ) *** " << endl;
00299
00300 //unused: Int_t tmpint;
00301 const char* tmpchar = 0;
00302
00303 if(r.Get("ListIn",tmpchar)){ fListIn = tmpchar; }
00304 if(r.Get("ListOut",tmpchar)){ fListOut = tmpchar; }
00305
00306 AlgFactory &af = AlgFactory::GetInstance();
00307
00308 AlgHandle ah_list = af.GetAlgHandle("AlgFarDetShieldPlankList","default");
00309 AlgConfig &ac_list = ah_list.GetAlgConfig();
00310 ac_list.UnLockValues();
00311 ac_list.LockValues();
00312
00313 AlgHandle ah_plank = af.GetAlgHandle("AlgFarDetShieldPlank","default");
00314 AlgConfig &ac_plank = ah_plank.GetAlgConfig();
00315 ac_plank.UnLockValues();
00316 ac_plank.LockValues();
00317
00318 return;
00319 }
|
|
|
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 268 of file FarDetShieldPlankListModule.cxx. References fListIn, fListOut, AlgHandle::GetAlgConfig(), AlgFactory::GetAlgHandle(), AlgFactory::GetInstance(), Registry::LockValues(), MSG, AlgFactory::Register(), Registry::Set(), and Registry::UnLockValues(). 00269 {
00270 MSG("FarDetShieldPlank",Msg::kDebug) << " *** FarDetShieldPlankListModule::DefaultConfig( ) *** " << endl;
00271
00272 static Registry r;
00273 r.SetName("FarDetShieldPlankListModule.config.default");
00274 r.UnLockValues();
00275 r.Set("ListIn",fListIn.Data());
00276 r.Set("ListOut",fListOut.Data());
00277 r.LockValues();
00278
00279 AlgFactory &af = AlgFactory::GetInstance();
00280
00281 af.Register("AlgFarDetShieldPlankList","default");
00282 AlgHandle ah_list = af.GetAlgHandle("AlgFarDetShieldPlankList","default");
00283 AlgConfig &ac_list = ah_list.GetAlgConfig();
00284 ac_list.UnLockValues();
00285 ac_list.LockValues();
00286
00287 af.Register("AlgFarDetShieldPlank","default");
00288 AlgHandle ah_plank = af.GetAlgHandle("AlgFarDetShieldPlank","default");
00289 AlgConfig &ac_plank = ah_plank.GetAlgConfig();
00290 ac_plank.UnLockValues();
00291 ac_plank.LockValues();
00292
00293 return r;
00294 }
|
|
|
Implement for notification of end of job Reimplemented from JobCModule. Definition at line 321 of file FarDetShieldPlankListModule.cxx. References fShieldFile, fShieldTree, and MSG. 00322 {
00323 if(fShieldFile){
00324 MSG("FarDetShieldPlank",Msg::kDebug) << " *** saving SHIELD INFO to file *** " << endl;
00325 TDirectory* tmpd = gDirectory;
00326 fShieldFile->cd();
00327 fShieldTree->Write();
00328 fShieldFile->Close();
00329 gDirectory = tmpd;
00330 MSG("FarDetShieldPlank",Msg::kDebug) << " ... data saved to file *** " << endl;
00331 }
00332 }
|
|
|
Implement this for read-write access to the MomNavigator Reimplemented from JobCModule. Definition at line 56 of file FarDetShieldPlankListModule.cxx. References CandRecord::FindCandHandle(), fListIn, fListOut, AlgFactory::GetAlgHandle(), MomNavigator::GetFragment(), AlgFactory::GetInstance(), FarDetShieldPlankList::MakeCandidate(), MSG, CandRecord::SecureCandHandle(), CandContext::SetCandRecord(), CandContext::SetDataIn(), JobCResult::SetFailed(), CandHandle::SetName(), and CandHandle::SetTitle(). 00057 {
00058 MSG("FarDetShieldPlank",Msg::kDebug) << " *** FarDetShieldPlankListModule::Reco(...) *** " << endl;
00059 JobCResult result(JobCResult::kPassed);
00060
00061 CandRecord *candrec = dynamic_cast<CandRecord*>(mom->GetFragment("CandRecord", "PrimaryCandidateRecord"));
00062 if (candrec == 0) {
00063 MSG("FarDetShieldPlank", Msg::kWarning) << " *** FAILED TO FIND CANDRECORD *** " << endl;
00064 return result.SetFailed();
00065 }
00066
00067 CandDigitListHandle *cdlh = dynamic_cast<CandDigitListHandle*>(candrec->FindCandHandle("CandDigitListHandle", fListIn.Data()));
00068 if(cdlh){
00069 MSG("FarDetShieldPlank", Msg::kDebug) << " *** MAKE HANDLE *** " << endl;
00070 TObjArray* myarray = new TObjArray();
00071 myarray->Add(cdlh);
00072 AlgFactory &af = AlgFactory::GetInstance();
00073 AlgHandle ah = af.GetAlgHandle("AlgFarDetShieldPlankList","default");
00074
00075 CandContext cx(this, mom);
00076 cx.SetCandRecord(candrec);
00077 cx.SetDataIn(myarray);
00078
00079 FarDetShieldPlankListHandle csslh = FarDetShieldPlankList::MakeCandidate(ah, cx);
00080 csslh.SetName(fListOut.Data());
00081 csslh.SetTitle(TString("Created by FarDetShieldPlankListModule"));
00082
00083 candrec->SecureCandHandle(csslh);
00084
00085 delete myarray;
00086 }
00087
00088 return result;
00089 }
|
|
|
Definition at line 41 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 53 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 31 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 24 of file FarDetShieldPlankListModule.h. Referenced by Config(), DefaultConfig(), and Reco(). |
|
|
Definition at line 25 of file FarDetShieldPlankListModule.h. Referenced by Config(), DefaultConfig(), and Reco(). |
|
|
Definition at line 65 of file FarDetShieldPlankListModule.h. |
|
|
Definition at line 66 of file FarDetShieldPlankListModule.h. |
|
|
Definition at line 36 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 37 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 35 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 34 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 50 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 62 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 49 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 61 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 27 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 32 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 28 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 33 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 48 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 60 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 30 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 29 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 47 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 59 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 44 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 56 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 46 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 58 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 45 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 57 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 42 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 54 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 43 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 55 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 38 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 39 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 51 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
|
|
Definition at line 63 of file FarDetShieldPlankListModule.h. Referenced by Ana(). |
1.3.9.1