CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/Fireworks/Core/src/FWMagField.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 
00003 #include "TH1F.h"
00004 #include "Fireworks/Core/interface/FWMagField.h"
00005 #include "Fireworks/Core/interface/fwLog.h"
00006 #include "DataFormats/FWLite/interface/Handle.h"
00007 #include "DataFormats/Common/interface/ConditionsInEdm.h"
00008 #include "DataFormats/Scalers/interface/DcsStatus.h"
00009 #include "DataFormats/FWLite/interface/Event.h"
00010 
00011 FWMagField::FWMagField() :
00012    TEveMagField(),
00013 
00014    m_source(kNone),
00015    m_userField(-1),
00016    m_eventField(-1),
00017 
00018    m_reverse(true),
00019    m_simpleModel(false),
00020 
00021    m_guessValHist(0),
00022    m_numberOfFieldIsOnEstimates(0),
00023    m_numberOfFieldEstimates(0),
00024    m_updateFieldEstimate(true),
00025    m_guessedField(0)
00026 {
00027    m_guessValHist = new TH1F("FieldEstimations", "Field estimations from tracks and muons",
00028                              200, -4.5, 4.5);
00029    m_guessValHist->SetDirectory(0);
00030 }
00031 
00032 FWMagField::~FWMagField()
00033 {
00034    delete m_guessValHist;
00035 }
00036 
00037 //______________________________________________________________________________
00038 
00039 TEveVector
00040 FWMagField::GetField(Float_t x, Float_t y, Float_t z) const
00041 {
00042    // Virtual method of TEveMagField class.
00043 
00044    Float_t R = sqrt(x*x+y*y);
00045    Float_t field = m_reverse ? -GetFieldMag() : GetFieldMag();
00046 
00047    //barrel
00048    if ( TMath::Abs(z)<724 ){
00049     
00050       //inside solenoid
00051       if ( R < 300) return TEveVector(0,0,field);
00052       // outside solinoid
00053       if ( m_simpleModel ||
00054            ( R>461.0 && R<490.5 ) ||
00055            ( R>534.5 && R<597.5 ) ||
00056            ( R>637.0 && R<700.0 ) )
00057          return TEveVector(0,0,-field/3.8*1.2);
00058   
00059    } else {
00060       // endcaps
00061       if (m_simpleModel){
00062          if ( R < 50 ) return TEveVector(0,0,field);
00063          if ( z > 0 )
00064             return TEveVector(x/R*field/3.8*2.0, y/R*field/3.8*2.0, 0);
00065          else
00066             return TEveVector(-x/R*field/3.8*2.0, -y/R*field/3.8*2.0, 0);
00067       }
00068       // proper model
00069       if ( ( ( TMath::Abs(z)>724 ) && ( TMath::Abs(z)<786 ) ) ||
00070            ( ( TMath::Abs(z)>850 ) && ( TMath::Abs(z)<910 ) ) ||
00071            ( ( TMath::Abs(z)>975 ) && ( TMath::Abs(z)<1003 ) ) )
00072       {
00073          if ( z > 0 )
00074             return TEveVector(x/R*field/3.8*2.0, y/R*field/3.8*2.0, 0);
00075          else
00076             return TEveVector(-x/R*field/3.8*2.0, -y/R*field/3.8*2.0, 0);
00077       }
00078    }
00079    return TEveVector(0,0,0);
00080 }
00081 
00082 //______________________________________________________________________________
00083 
00084 Float_t
00085 FWMagField::GetFieldMag() const
00086 {
00087    float res;
00088    switch ( m_source )
00089    {
00090       case kEvent:
00091       {
00092          res = m_eventField;  
00093          break;
00094       }
00095       case kUser:
00096       {
00097          res = m_userField;   
00098          break;
00099       }
00100       default:
00101       {
00102          if ( m_updateFieldEstimate )
00103          {
00104             if ( m_guessValHist->GetEntries() > 2  && m_guessValHist->GetRMS()  < 0.5 )
00105             {
00106                m_guessedField = m_guessValHist->GetMean();
00107             
00108 
00109                // std::cout << "FWMagField::GetFieldMag(), get average "
00110                //  << m_guessValHist->GetMean() << " guessed value: RMS= "<< m_guessValHist->GetRMS()
00111                //  <<" samples "<< m_guessValHist->GetEntries() << std::endl;
00112             
00113             }
00114             else if ( m_numberOfFieldIsOnEstimates > m_numberOfFieldEstimates/2 || m_numberOfFieldEstimates == 0 )
00115             {
00116                m_guessedField = 3.8;
00117                // fwLog(fwlog::kDebug) << "FWMagField::GetFieldMag() get default field, number estimates "
00118                //  << m_numberOfFieldEstimates << " number fields is on  m_numberOfFieldIsOnEstimates" <<std::endl;
00119             }
00120             else
00121             {
00122                m_guessedField = 0;
00123                //  fwLog(fwlog::kDebug) << "Update field estimate, guess field is OFF." <<std::endl;
00124             }
00125             m_updateFieldEstimate  = false;
00126          }
00127 
00128          res = m_guessedField;
00129       }
00130    }
00131 
00132    return res;
00133 }
00134 
00135 Float_t
00136 FWMagField::GetMaxFieldMag() const
00137 {
00138    // Runge-Kutta stepper does not like this to be zero.
00139    // Will be fixed in root.
00140    // The return value should definitley not be negative -- so Abs
00141    // should stay!
00142 
00143    return TMath::Max(TMath::Abs(GetFieldMag()), 0.01f);
00144 }
00145 
00146 //______________________________________________________________________________
00147 
00148 void FWMagField::guessFieldIsOn(bool isOn) const
00149 {
00150    if ( isOn ) ++m_numberOfFieldIsOnEstimates;
00151    ++m_numberOfFieldEstimates;
00152    m_updateFieldEstimate  = true;
00153 }
00154 
00155 void FWMagField::guessField(float val) const
00156 {
00157    // fwLog(filedDebug) <<  "FWMagField::guessField "<< val << std::endl;
00158    m_guessValHist->Fill(val);
00159    m_updateFieldEstimate = true; 
00160 }
00161 
00162 void FWMagField::resetFieldEstimate() const
00163 {
00164    m_guessValHist->Reset();
00165    m_guessValHist->SetAxisRange(-4, 4);
00166    m_numberOfFieldIsOnEstimates = 0;
00167    m_numberOfFieldEstimates = 0;
00168    m_updateFieldEstimate = true;   
00169 }
00170 
00171 //______________________________________________________________________________
00172 void FWMagField::checkFieldInfo(const edm::EventBase* event)
00173 {
00174    const static float  currentToField = 3.8/18160;
00175    bool available = false;
00176    try
00177    {
00178       edm::InputTag conditionsTag("conditionsInEdm");
00179       edm::Handle<edm::ConditionsInRunBlock> runCond;
00180       // FIXME: ugly hack to avoid exposing an fwlite::Event. Need to ask
00181       //        Chris / Ianna how to get mag field from an EventBase.
00182       const fwlite::Event *fwEvent = dynamic_cast<const fwlite::Event*>(event);
00183       if (!fwEvent)
00184          return;
00185       fwEvent->getRun().getByLabel(conditionsTag, runCond);
00186       
00187       if( runCond.isValid())
00188       {
00189          available = true;
00190          m_eventField = currentToField * runCond->BAvgCurrent;
00191          fwLog( fwlog::kDebug ) << "Magnetic field info found in ConditionsInEdm branch : "<< m_eventField << std::endl;
00192       }
00193       else
00194       {
00195          edm::InputTag dcsTag("scalersRawToDigi");
00196          edm::Handle< std::vector<DcsStatus> > dcsStatus;
00197          event->getByLabel(dcsTag, dcsStatus);
00198          
00199          if (dcsStatus.isValid() && dcsStatus->size())
00200          {
00201             float sum = 0;
00202             for (std::vector<DcsStatus>::const_iterator i = dcsStatus->begin(); i <  dcsStatus->end(); ++i)
00203                sum += (*i).magnetCurrent();
00204 
00205             available = true;
00206             m_eventField = currentToField * sum/dcsStatus->size();
00207             fwLog( fwlog::kDebug) << "Magnetic field info found in DcsStatus branch: " << m_eventField << std::endl;
00208          }
00209 
00210       }
00211    }
00212    catch (cms::Exception&)
00213    {
00214       fwLog( fwlog::kDebug ) << "Cought exception in FWMagField::checkFieldInfo\n";
00215    }
00216 
00217    if (!available)
00218    {
00219       m_source = kNone;
00220       fwLog( fwlog::kDebug ) << "No magnetic field info available in Event\n";
00221    }
00222 }
00223