![]() |
![]() |
00001 00002 // $Id: LumiSummary.cc,v 1.25 2012/03/06 11:39:18 xiezhen Exp $ 00003 00004 #include "DataFormats/Luminosity/interface/LumiSummary.h" 00005 00006 #include <iomanip> 00007 #include <ostream> 00008 #include <iostream> 00009 float 00010 LumiSummary::avgInsDelLumi()const{ 00011 if(lumiversion_=="v2"){ 00012 return avginsdellumi_*1000.0; 00013 } 00014 return avginsdellumi_; 00015 } 00016 float 00017 LumiSummary::intgDelLumi()const{ 00018 return this->avgInsDelLumi()*float(this->lumiSectionLength()); 00019 } 00020 float 00021 LumiSummary::avgInsDelLumiErr()const{ 00022 return avginsdellumierr_; 00023 } 00024 float 00025 LumiSummary::intgRecLumi()const{ 00026 return this->avgInsRecLumi()*float(this->lumiSectionLength()); 00027 } 00028 short 00029 LumiSummary::lumiSecQual()const { 00030 return lumisecqual_; 00031 } 00032 unsigned long long 00033 LumiSummary::deadcount() const{ 00034 return deadcount_; 00035 } 00036 unsigned long long 00037 LumiSummary::bitzerocount() const{ 00038 return bitzerocount_; 00039 } 00040 float 00041 LumiSummary::deadFrac() const { 00042 //definition: deadcount/bitzerocount 00043 //if no trigger data, return deadfraction 1.0,mask out this LS 00044 //if bitzerocount=0, return -1.0 meaning no beam 00045 if (l1data_.size()==0) return 1.0; 00046 if (bitzerocount_==0) return -1.0; 00047 return float(deadcount_)/float(bitzerocount_); 00048 } 00049 float 00050 LumiSummary::liveFrac() const { 00051 //1-deadfraction 00052 //else if deadfraction<0 meaning no beam, live fraction=0 00053 // 00054 if (deadFrac()<0) return 0; 00055 return 1-deadFrac(); 00056 } 00057 float 00058 LumiSummary::lumiSectionLength() const { 00059 //numorbits*3564*24.95e-09 00060 return numorbit_*3564.0*24.95e-9; 00061 } 00062 unsigned int 00063 LumiSummary::lsNumber() const{ 00064 return lsnumber_; 00065 } 00066 unsigned int 00067 LumiSummary::startOrbit() const{ 00068 return startorbit_; 00069 } 00070 unsigned int 00071 LumiSummary::numOrbit() const{ 00072 return numorbit_; 00073 } 00074 bool 00075 LumiSummary::isValid() const { 00076 return (lumiversion_!="-1"); 00077 } 00078 LumiSummary::L1 00079 LumiSummary::l1info(unsigned int idx)const{ 00080 return l1data_.at(idx); 00081 } 00082 LumiSummary::HLT 00083 LumiSummary::hltinfo(unsigned int idx) const { 00084 return hltdata_.at(idx); 00085 } 00086 size_t 00087 LumiSummary::nTriggerLine()const{ 00088 return l1data_.size(); 00089 } 00090 size_t 00091 LumiSummary::nHLTPath()const{ 00092 return hltdata_.size(); 00093 } 00094 float 00095 LumiSummary::avgInsRecLumi() const { 00096 return this->avgInsDelLumi() * liveFrac(); 00097 } 00098 float 00099 LumiSummary::avgInsRecLumiErr() const { 00100 return avginsdellumierr_ * liveFrac(); 00101 } 00102 bool 00103 LumiSummary::isProductEqual(LumiSummary const& next) const { 00104 return (avginsdellumi_ == next.avginsdellumi_ && 00105 avginsdellumierr_ == next.avginsdellumierr_ && 00106 lumisecqual_ == next.lumisecqual_ && 00107 deadcount_ == next.deadcount_ && 00108 lsnumber_ == next.lsnumber_ && 00109 startorbit_== next.startorbit_ && 00110 numorbit_==next.numorbit_&& 00111 l1data_.size() == next.l1data_.size() && 00112 hltdata_.size() == next.hltdata_.size() && 00113 lumiversion_ == next.lumiversion_ ); 00114 } 00115 std::string 00116 LumiSummary::lumiVersion()const{ 00117 return lumiversion_; 00118 } 00119 void 00120 LumiSummary::setLumiVersion(const std::string& lumiversion){ 00121 lumiversion_=lumiversion; 00122 } 00123 void 00124 LumiSummary::setLumiData(float instlumi,float instlumierr,short lumiquality){ 00125 avginsdellumi_=instlumi; 00126 avginsdellumierr_=instlumierr; 00127 lumisecqual_=lumiquality; 00128 } 00129 void 00130 LumiSummary::setDeadCount(unsigned long long deadcount){ 00131 deadcount_=deadcount; 00132 } 00133 void 00134 LumiSummary::setBitZeroCount(unsigned long long bitzerocount){ 00135 bitzerocount_=bitzerocount; 00136 } 00137 void 00138 LumiSummary::setlsnumber(unsigned int lsnumber){ 00139 lsnumber_=lsnumber; 00140 } 00141 void 00142 LumiSummary::setOrbitData(unsigned int startorbit,unsigned int numorbit){ 00143 startorbit_=startorbit; 00144 numorbit_=numorbit; 00145 } 00146 void 00147 LumiSummary::swapL1Data(std::vector<L1>& l1data){ 00148 l1data_.swap(l1data); 00149 } 00150 void 00151 LumiSummary::swapHLTData(std::vector<HLT>& hltdata){ 00152 hltdata_.swap(hltdata); 00153 } 00154 void 00155 LumiSummary::copyL1Data(const std::vector<L1>& l1data){ 00156 l1data_.assign(l1data.begin(),l1data.end()); 00157 } 00158 void 00159 LumiSummary::copyHLTData(const std::vector<HLT>& hltdata){ 00160 hltdata_.assign(hltdata.begin(),hltdata.end()); 00161 } 00162 std::ostream& operator<<(std::ostream& s, const LumiSummary& lumiSummary) { 00163 s << "\nDumping LumiSummary\n\n"; 00164 if(!lumiSummary.isValid()){ 00165 s << " === Invalid Lumi values === \n"; 00166 } 00167 s << " lumiVersion = " << lumiSummary.lumiVersion() << "\n"; 00168 s << " avgInsDelLumi = " << lumiSummary.avgInsDelLumi() << "\n"; 00169 s << " avgIntgDelLumi = " << lumiSummary.intgDelLumi() <<"\n"; 00170 s << " avgInsDelLumiErr = " << lumiSummary.avgInsDelLumiErr() << "\n"; 00171 s << " lumiSecQual = " << lumiSummary.lumiSecQual() << "\n"; 00172 s << " deadCount = " << lumiSummary.deadcount() << "\n"; 00173 s << " bitZeroCount = " << lumiSummary.bitzerocount() << "\n"; 00174 s << " deadFrac = " << (float)lumiSummary.deadFrac() << "\n"; 00175 s << " liveFrac = " << (float)lumiSummary.liveFrac() << "\n"; 00176 s << " lsNumber = " << lumiSummary.lsNumber() << "\n"; 00177 s << " startOrbit = " << lumiSummary.startOrbit() <<"\n"; 00178 s << " numOrbit = " << lumiSummary.numOrbit() <<"\n"; 00179 s << " avgInsRecLumi = " << lumiSummary.avgInsRecLumi() << "\n"; 00180 s << " avgInsRecLumiErr = " << lumiSummary.avgInsRecLumiErr() << "\n\n"; 00181 s << std::setw(15) << "l1nameidx"; 00182 s << std::setw(15) << "l1prescale"; 00183 s << "\n"; 00184 size_t nTriggers=lumiSummary.nTriggerLine(); 00185 size_t nHLTPath=lumiSummary.nHLTPath(); 00186 for(unsigned int i = 0; i < nTriggers; ++i) { 00187 s << std::setw(15); 00188 s << lumiSummary.l1info(i).triggernameidx; 00189 s << std::setw(15); 00190 s << lumiSummary.l1info(i).prescale; 00191 s<<"\n"; 00192 } 00193 s << std::setw(15) << "hltpathidx"; 00194 s << std::setw(15) << "hltprescale"; 00195 s << "\n"; 00196 for(unsigned int i = 0; i < nHLTPath; ++i) { 00197 s << std::setw(15); 00198 s << lumiSummary.hltinfo(i).pathnameidx; 00199 s << std::setw(15); 00200 s << lumiSummary.hltinfo(i).prescale; 00201 s << "\n"; 00202 } 00203 return s << "\n"; 00204 }