#include <FWMagField.h>
Public Types | |
enum | ESource { kNone, kEvent, kUser } |
Public Member Functions | |
void | checkFieldInfo (const edm::EventBase *) |
FWMagField () | |
virtual TEveVector | GetField (Float_t x, Float_t y, Float_t z) const |
virtual Float_t | GetMaxFieldMag () const |
ESource | getSource () const |
float | getUserField () const |
void | guessField (float estimate) const |
void | guessFieldIsOn (bool guess) const |
bool | isReverse () const |
bool | isSimpleModel () const |
void | resetFieldEstimate () const |
void | setReverseState (bool state) |
void | setSimpleModel (bool simpleModel) |
void | setSource (ESource x) |
void | setUserField (float b) |
virtual | ~FWMagField () |
Private Member Functions | |
FWMagField (const FWMagField &) | |
Float_t | GetFieldMag () const |
const FWMagField & | operator= (const FWMagField &) |
Private Attributes | |
float | m_eventField |
float | m_guessedField |
TH1F * | m_guessValHist |
int | m_numberOfFieldEstimates |
int | m_numberOfFieldIsOnEstimates |
bool | m_reverse |
bool | m_simpleModel |
ESource | m_source |
bool | m_updateFieldEstimate |
float | m_userField |
Definition at line 16 of file FWMagField.h.
enum FWMagField::ESource |
Definition at line 21 of file FWMagField.h.
FWMagField::FWMagField | ( | ) |
Definition at line 11 of file FWMagField.cc.
References m_guessValHist.
: TEveMagField(), m_source(kNone), m_userField(-1), m_eventField(-1), m_reverse(true), m_simpleModel(false), m_guessValHist(0), m_numberOfFieldIsOnEstimates(0), m_numberOfFieldEstimates(0), m_updateFieldEstimate(true), m_guessedField(0) { m_guessValHist = new TH1F("FieldEstimations", "Field estimations from tracks and muons", 200, -4.5, 4.5); m_guessValHist->SetDirectory(0); }
FWMagField::~FWMagField | ( | ) | [virtual] |
FWMagField::FWMagField | ( | const FWMagField & | ) | [private] |
void FWMagField::checkFieldInfo | ( | const edm::EventBase * | event | ) |
Definition at line 172 of file FWMagField.cc.
References event(), fwLog, fwlite::Run::getByLabel(), fwlite::Event::getRun(), i, edm::HandleBase::isValid(), fwlog::kDebug, kNone, m_eventField, and m_source.
Referenced by CmsShowMainBase::draw().
{ const static float currentToField = 3.8/18160; bool available = false; try { edm::InputTag conditionsTag("conditionsInEdm"); edm::Handle<edm::ConditionsInRunBlock> runCond; // FIXME: ugly hack to avoid exposing an fwlite::Event. Need to ask // Chris / Ianna how to get mag field from an EventBase. const fwlite::Event *fwEvent = dynamic_cast<const fwlite::Event*>(event); if (!fwEvent) return; fwEvent->getRun().getByLabel(conditionsTag, runCond); if( runCond.isValid()) { available = true; m_eventField = currentToField * runCond->BAvgCurrent; fwLog( fwlog::kDebug ) << "Magnetic field info found in ConditionsInEdm branch : "<< m_eventField << std::endl; } else { edm::InputTag dcsTag("scalersRawToDigi"); edm::Handle< std::vector<DcsStatus> > dcsStatus; event->getByLabel(dcsTag, dcsStatus); if (dcsStatus.isValid() && dcsStatus->size()) { float sum = 0; for (std::vector<DcsStatus>::const_iterator i = dcsStatus->begin(); i < dcsStatus->end(); ++i) sum += (*i).magnetCurrent(); available = true; m_eventField = currentToField * sum/dcsStatus->size(); fwLog( fwlog::kDebug) << "Magnetic field info found in DcsStatus branch: " << m_eventField << std::endl; } } } catch (cms::Exception&) { fwLog( fwlog::kDebug ) << "Cought exception in FWMagField::checkFieldInfo\n"; } if (!available) { m_source = kNone; fwLog( fwlog::kDebug ) << "No magnetic field info available in Event\n"; } }
TEveVector FWMagField::GetField | ( | Float_t | x, |
Float_t | y, | ||
Float_t | z | ||
) | const [virtual] |
Definition at line 40 of file FWMagField.cc.
References GetFieldMag(), m_reverse, m_simpleModel, dttmaxenums::R, and mathSSE::sqrt().
{ // Virtual method of TEveMagField class. Float_t R = sqrt(x*x+y*y); Float_t field = m_reverse ? -GetFieldMag() : GetFieldMag(); //barrel if ( TMath::Abs(z)<724 ){ //inside solenoid if ( R < 300) return TEveVector(0,0,field); // outside solinoid if ( m_simpleModel || ( R>461.0 && R<490.5 ) || ( R>534.5 && R<597.5 ) || ( R>637.0 && R<700.0 ) ) return TEveVector(0,0,-field/3.8*1.2); } else { // endcaps if (m_simpleModel){ if ( R < 50 ) return TEveVector(0,0,field); if ( z > 0 ) return TEveVector(x/R*field/3.8*2.0, y/R*field/3.8*2.0, 0); else return TEveVector(-x/R*field/3.8*2.0, -y/R*field/3.8*2.0, 0); } // proper model if ( ( ( TMath::Abs(z)>724 ) && ( TMath::Abs(z)<786 ) ) || ( ( TMath::Abs(z)>850 ) && ( TMath::Abs(z)<910 ) ) || ( ( TMath::Abs(z)>975 ) && ( TMath::Abs(z)<1003 ) ) ) { if ( z > 0 ) return TEveVector(x/R*field/3.8*2.0, y/R*field/3.8*2.0, 0); else return TEveVector(-x/R*field/3.8*2.0, -y/R*field/3.8*2.0, 0); } } return TEveVector(0,0,0); }
Float_t FWMagField::GetFieldMag | ( | ) | const [private] |
Definition at line 85 of file FWMagField.cc.
References kEvent, kUser, m_eventField, m_guessedField, m_guessValHist, m_numberOfFieldEstimates, m_numberOfFieldIsOnEstimates, m_source, m_updateFieldEstimate, and m_userField.
Referenced by GetField(), and GetMaxFieldMag().
{ float res; switch ( m_source ) { case kEvent: { res = m_eventField; break; } case kUser: { res = m_userField; break; } default: { if ( m_updateFieldEstimate ) { if ( m_guessValHist->GetEntries() > 2 && m_guessValHist->GetRMS() < 0.5 ) { m_guessedField = m_guessValHist->GetMean(); // std::cout << "FWMagField::GetFieldMag(), get average " // << m_guessValHist->GetMean() << " guessed value: RMS= "<< m_guessValHist->GetRMS() // <<" samples "<< m_guessValHist->GetEntries() << std::endl; } else if ( m_numberOfFieldIsOnEstimates > m_numberOfFieldEstimates/2 || m_numberOfFieldEstimates == 0 ) { m_guessedField = 3.8; // fwLog(fwlog::kDebug) << "FWMagField::GetFieldMag() get default field, number estimates " // << m_numberOfFieldEstimates << " number fields is on m_numberOfFieldIsOnEstimates" <<std::endl; } else { m_guessedField = 0; // fwLog(fwlog::kDebug) << "Update field estimate, guess field is OFF." <<std::endl; } m_updateFieldEstimate = false; } res = m_guessedField; } } return res; }
Float_t FWMagField::GetMaxFieldMag | ( | ) | const [virtual] |
Definition at line 136 of file FWMagField.cc.
References f, GetFieldMag(), and siStripFEDMonitor_P5_cff::Max.
{ // Runge-Kutta stepper does not like this to be zero. // Will be fixed in root. // The return value should definitley not be negative -- so Abs // should stay! return TMath::Max(TMath::Abs(GetFieldMag()), 0.01f); }
ESource FWMagField::getSource | ( | ) | const [inline] |
Definition at line 34 of file FWMagField.h.
References m_source.
Referenced by FWMuonBuilder::calculateField(), CmsShowMainBase::draw(), FWPFTrackUtils::setupLegoTrack(), and FWPFTrackUtils::setupTrack().
{ return m_source; }
float FWMagField::getUserField | ( | ) | const [inline] |
void FWMagField::guessField | ( | float | estimate | ) | const |
Definition at line 155 of file FWMagField.cc.
References m_guessValHist, and m_updateFieldEstimate.
Referenced by FWTrackProxyBuilder::build(), FWMuonBuilder::calculateField(), FWPFTrackUtils::setupLegoTrack(), and FWPFTrackUtils::setupTrack().
{ // fwLog(filedDebug) << "FWMagField::guessField "<< val << std::endl; m_guessValHist->Fill(val); m_updateFieldEstimate = true; }
void FWMagField::guessFieldIsOn | ( | bool | guess | ) | const |
Definition at line 148 of file FWMagField.cc.
References m_numberOfFieldEstimates, m_numberOfFieldIsOnEstimates, and m_updateFieldEstimate.
Referenced by FWMuonBuilder::calculateField().
{ if ( isOn ) ++m_numberOfFieldIsOnEstimates; ++m_numberOfFieldEstimates; m_updateFieldEstimate = true; }
bool FWMagField::isReverse | ( | ) | const [inline] |
bool FWMagField::isSimpleModel | ( | ) | const [inline] |
const FWMagField& FWMagField::operator= | ( | const FWMagField & | ) | [private] |
void FWMagField::resetFieldEstimate | ( | ) | const |
Definition at line 162 of file FWMagField.cc.
References m_guessValHist, m_numberOfFieldEstimates, m_numberOfFieldIsOnEstimates, and m_updateFieldEstimate.
{ m_guessValHist->Reset(); m_guessValHist->SetAxisRange(-4, 4); m_numberOfFieldIsOnEstimates = 0; m_numberOfFieldEstimates = 0; m_updateFieldEstimate = true; }
void FWMagField::setReverseState | ( | bool | state | ) | [inline] |
Definition at line 37 of file FWMagField.h.
References m_reverse, and evf::utils::state.
void FWMagField::setSimpleModel | ( | bool | simpleModel | ) | [inline] |
Definition at line 39 of file FWMagField.h.
References m_simpleModel.
{ m_simpleModel = simpleModel; }
void FWMagField::setSource | ( | ESource | x | ) | [inline] |
void FWMagField::setUserField | ( | float | b | ) | [inline] |
float FWMagField::m_eventField [private] |
Definition at line 55 of file FWMagField.h.
Referenced by checkFieldInfo(), and GetFieldMag().
float FWMagField::m_guessedField [mutable, private] |
Definition at line 65 of file FWMagField.h.
Referenced by GetFieldMag().
TH1F* FWMagField::m_guessValHist [mutable, private] |
Definition at line 61 of file FWMagField.h.
Referenced by FWMagField(), GetFieldMag(), guessField(), resetFieldEstimate(), and ~FWMagField().
int FWMagField::m_numberOfFieldEstimates [mutable, private] |
Definition at line 63 of file FWMagField.h.
Referenced by GetFieldMag(), guessFieldIsOn(), and resetFieldEstimate().
int FWMagField::m_numberOfFieldIsOnEstimates [mutable, private] |
Definition at line 62 of file FWMagField.h.
Referenced by GetFieldMag(), guessFieldIsOn(), and resetFieldEstimate().
bool FWMagField::m_reverse [private] |
Definition at line 57 of file FWMagField.h.
Referenced by GetField(), isReverse(), and setReverseState().
bool FWMagField::m_simpleModel [private] |
Definition at line 58 of file FWMagField.h.
Referenced by GetField(), isSimpleModel(), and setSimpleModel().
ESource FWMagField::m_source [private] |
Definition at line 53 of file FWMagField.h.
Referenced by checkFieldInfo(), GetFieldMag(), getSource(), and setSource().
bool FWMagField::m_updateFieldEstimate [mutable, private] |
Definition at line 64 of file FWMagField.h.
Referenced by GetFieldMag(), guessField(), guessFieldIsOn(), and resetFieldEstimate().
float FWMagField::m_userField [private] |
Definition at line 54 of file FWMagField.h.
Referenced by GetFieldMag(), getUserField(), and setUserField().