00001
00002
00003
00004
00005
00006
00007 class StJets;
00008
00009 #include "StEmcUtil/geometry/StEmcGeom.h"
00010 #include "StEEmcUtil/EEmcGeom/EEmcGeomSimple.h"
00011 #include "StEventTypes.h"
00012 #include "StMuDSTMaker/COMMON/StMuTypes.hh"
00013 #include "StMuTrackFourVec.h"
00014
00015 #include "StMuTrackEmu.h"
00016 #include "StMuTowerEmu.h"
00017 #include "StMcTrackEmu.h"
00018 #include "StSpinPool/StJetEvent/StJetEventTypes.h"
00019
00020 #include "StJetFinder/StProtoJet.h"
00021
00022 #include "StFourPMaker.h"
00023 #include "StBET4pMaker.h"
00024
00025 #include "TFile.h"
00026 #include "TTree.h"
00027
00028 #include "StjeJetEventTreeWriter.h"
00029
00030 ClassImp(StjeJetEventTreeWriter);
00031
00032 StjeJetEventTreeWriter::StjeJetEventTreeWriter(const char* outFileName)
00033 : _OutFileName(outFileName)
00034 , _jetTree(0)
00035 , _outFile(0)
00036 {
00037 }
00038
00039 void StjeJetEventTreeWriter::addJetFinder(StFourPMaker* fourPMaker, const vector<const AbstractFourVec*>* particleList, list<StProtoJet>* protoJetList, const char* name, StJets*)
00040 {
00041 AnalyzerCtl anaCtl;
00042
00043 anaCtl._branchName = name;
00044 anaCtl._fourPMaker = fourPMaker;
00045 anaCtl._protoJetList = protoJetList;
00046 anaCtl._jetEvent = new StJetEvent;
00047
00048 _analyzerCtlList.push_back(anaCtl);
00049 }
00050
00051 void StjeJetEventTreeWriter::Init()
00052 {
00053 _outFile = new TFile(_OutFileName, "recreate");
00054 _jetTree = new TTree("jet", "jetTree");
00055
00056 for (vector<AnalyzerCtl>::iterator iAnalyzer = _analyzerCtlList.begin(); iAnalyzer != _analyzerCtlList.end(); ++iAnalyzer)
00057 _jetTree->Branch(iAnalyzer->_branchName.c_str(), "StJetEvent", &iAnalyzer->_jetEvent);
00058
00059 _jetTree->BranchRef();
00060 }
00061
00062 void StjeJetEventTreeWriter::Finish()
00063 {
00064 _outFile->Write();
00065 _outFile->Close();
00066 delete _outFile;
00067 _outFile = 0;
00068 }
00069
00070 void StjeJetEventTreeWriter::fillJetTreeHeader(int iAnalyzer)
00071 {
00072 StJetEvent* jetEvent = _analyzerCtlList[iAnalyzer]._jetEvent;
00073 StFourPMaker* fourPMaker = _analyzerCtlList[iAnalyzer]._fourPMaker;
00074 jetEvent->Clear();
00075 jetEvent->mRunId = fourPMaker->GetRunNumber();
00076 jetEvent->mEventId = fourPMaker->GetEventNumber();
00077 jetEvent->mDatime = fourPMaker->GetDateTime();
00078 }
00079
00080 void StjeJetEventTreeWriter::fillJetTree(int iAnalyzer, int iVertex)
00081 {
00082 AnalyzerCtl& analyzerCtl = _analyzerCtlList[iAnalyzer];
00083 fillJetTreeForOneVertex(analyzerCtl._jetEvent,analyzerCtl._protoJetList,analyzerCtl._fourPMaker,iVertex);
00084 }
00085
00086 void StjeJetEventTreeWriter::fillJetTreeForOneVertex(StJetEvent* jetEvent, list<StProtoJet>* protoJetList, StFourPMaker* fourPMaker, int iVertex)
00087 {
00088 StJetVertex* jetVertex = jetEvent->newVertex();
00089 copyVertex(fourPMaker->getVertexNodes()[iVertex].vertex,jetVertex);
00090 for (list<StProtoJet>::iterator protojet = protoJetList->begin(); protojet != protoJetList->end(); ++protojet)
00091 jetVertex->addJet(fillJet(jetEvent,jetVertex,*protojet));
00092 }
00093
00094 StJetCandidate* StjeJetEventTreeWriter::fillJet(StJetEvent* jetEvent, StJetVertex* jetVertex, StProtoJet& protojet)
00095 {
00096 StJetCandidate* jet = jetEvent->newJet(jetVertex->position(),TLorentzVector(protojet.px(),protojet.py(),protojet.pz(),protojet.e()));
00097
00098
00099 const StProtoJet::FourVecList& particleList = protojet.list();
00100 for (StProtoJet::FourVecList::const_iterator iParticle = particleList.begin(); iParticle != particleList.end(); ++iParticle) {
00101 const StMuTrackFourVec* particle = dynamic_cast<const StMuTrackFourVec*>(*iParticle);
00102
00103 if (StMuTrackEmu* t = particle->track()) {
00104 StJetTrack* track = jetEvent->newTrack();
00105 track->mId = t->id();
00106 track->mDetectorId = t->detectorId();
00107 track->mFlag = t->flag();
00108 track->mCharge = t->charge();
00109 track->mNHits = t->nHits();
00110 track->mNHitsFit = t->nHitsFit();
00111 track->mNHitsPoss = t->nHitsPoss();
00112 track->mNHitsDedx = t->nHitsDedx();
00113 track->mDedx = t->dEdx();
00114 track->mBeta = t->beta();
00115 track->mFirstPoint = t->firstPoint();
00116 track->mLastPoint = t->lastPoint();
00117 track->mExitTowerId = t->exitTowerId();
00118 track->mExitDetectorId = t->exitDetectorId();
00119 track->mExitPoint.SetPtEtaPhi(t->bemcRadius(),t->etaext(),t->phiext());
00120 track->mDca.SetXYZ(t->dcaX(),t->dcaY(),t->dcaZ());
00121 track->mDcaD = t->dcaD();
00122 track->mChi2 = t->chi2();
00123 track->mChi2Prob = t->chi2prob();
00124 TVector3 mom(t->px(),t->py(),t->pz());
00125 track->mPt = mom.Pt();
00126 track->mEta = mom.Eta();
00127 track->mPhi = mom.Phi();
00128 jet->addTrack(track)->setJet(jet);
00129 }
00130
00131 if (StMuTowerEmu* t = particle->tower()) {
00132 StJetTower* tower = jetEvent->newTower();
00133 tower->mId = t->id();
00134 tower->mDetectorId = t->detectorId();
00135 tower->mAdc = t->adc();
00136 tower->mPedestal = t->pedestal();
00137 tower->mRms = t->rms();
00138 tower->mStatus = t->status();
00139
00140
00141
00142
00143 TVector3 mom(t->px(),t->py(),t->pz());
00144 float energy = mom.Mag();
00145
00146 switch (t->detectorId()) {
00147 case kBarrelEmcTowerId:
00148 mom.SetPtEtaPhi(StEmcGeom::instance("bemc")->Radius(),mom.Eta(),mom.Phi());
00149 break;
00150 case kEndcapEmcTowerId:
00151 mom.SetMag(EEmcGeomSimple::Instance().getZMean()/mom.Unit().z());
00152 break;
00153 }
00154
00155 mom -= jetVertex->position();
00156 mom.SetMag(energy);
00157
00158 tower->mPt = mom.Pt();
00159 tower->mEta = mom.Eta();
00160 tower->mPhi = mom.Phi();
00161 jet->addTower(tower)->setJet(jet);
00162 }
00163
00164 if (StMcTrackEmu* t = particle->mctrack()) {
00165 StJetParticle* part = jetEvent->newParticle();
00166 part->mId = t->id();
00167 part->mPt = t->pt();
00168 part->mEta = t->eta();
00169 part->mPhi = t->phi();
00170 part->mM = t->m();
00171 part->mE = t->e();
00172 part->mPdg = t->pdg();
00173 part->mStatus = t->status();
00174 jet->addParticle(part)->setJet(jet);
00175 }
00176 }
00177
00178 float sumTowerPt = jet->sumTowerPt();
00179 float sumTrackPt = jet->sumTrackPt();
00180
00181 jet->mRt = sumTowerPt/(sumTowerPt+sumTrackPt);
00182
00183 jet->setVertex(jetVertex);
00184
00185 return jet;
00186 }
00187
00188 void StjeJetEventTreeWriter::copyVertex(const StMuPrimaryVertex* muVertex, StJetVertex* jetVertex)
00189 {
00190 jetVertex->mPosition = muVertex->position().xyz();
00191 jetVertex->mPosError = muVertex->posError().xyz();
00192 jetVertex->mVertexFinderId = muVertex->vertexFinderId();
00193 jetVertex->mRanking = muVertex->ranking();
00194 jetVertex->mNTracksUsed = muVertex->nTracksUsed();
00195 jetVertex->mNBTOFMatch = muVertex->nBTOFMatch();
00196 jetVertex->mNCTBMatch = muVertex->nCTBMatch();
00197 jetVertex->mNBEMCMatch = muVertex->nBEMCMatch();
00198 jetVertex->mNEEMCMatch = muVertex->nEEMCMatch();
00199 jetVertex->mNCrossCentralMembrane = muVertex->nCrossCentralMembrane();
00200 jetVertex->mSumTrackPt = muVertex->sumTrackPt();
00201 jetVertex->mMeanDip = muVertex->meanDip();
00202 jetVertex->mChiSquared = muVertex->chiSquared();
00203 jetVertex->mRefMultPos = muVertex->refMultPos();
00204 jetVertex->mRefMultNeg = muVertex->refMultNeg();
00205 jetVertex->mRefMultFtpcWest = muVertex->refMultFtpcWest();
00206 jetVertex->mRefMultFtpcEast = muVertex->refMultFtpcEast();
00207 }