CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes

MaterialBudgetForward Class Reference

#include <MaterialBudgetForward.h>

Inheritance diagram for MaterialBudgetForward:
SimWatcher Observer< const BeginOfRun * > Observer< const BeginOfTrack * > Observer< const G4Step * > Observer< const EndOfTrack * >

List of all members.

Public Member Functions

 MaterialBudgetForward (const edm::ParameterSet &)
virtual ~MaterialBudgetForward ()

Private Member Functions

void book (const edm::ParameterSet &)
 MaterialBudgetForward (const MaterialBudgetForward &)
const MaterialBudgetForwardoperator= (const MaterialBudgetForward &)
bool stopAfter (const G4Step *)
void update (const G4Step *)
 This routine will be called when the appropriate signal arrives.
void update (const BeginOfRun *)
 This routine will be called when the appropriate signal arrives.
void update (const BeginOfTrack *)
 This routine will be called when the appropriate signal arrives.
void update (const EndOfTrack *)
 This routine will be called when the appropriate signal arrives.

Private Attributes

std::vector< double > boundaries
std::vector< int > constituents
std::vector< int > detLevels
std::vector< std::string > detNames
std::vector< std::string > detTypes
double eta
std::vector< double > etaRegions
std::vector< double > intLen
std::vector< G4LogicalVolume * > logVolumes
TProfile * me100 [maxSet]
TProfile * me200 [maxSet]
TProfile * me300 [maxSet]
TH1F * me400 [maxSet]
TProfile2D * me500 [maxSet]
TProfile2D * me600 [maxSet]
TProfile2D * me700 [maxSet]
TH2F * me800 [maxSet]
double phi
std::vector< double > radLen
std::vector< int > regionTypes
std::vector< int > stackOrder
std::vector< double > stepLen
double stepT

Static Private Attributes

static const int maxSet = 25

Detailed Description

Definition at line 24 of file MaterialBudgetForward.h.


Constructor & Destructor Documentation

MaterialBudgetForward::MaterialBudgetForward ( const edm::ParameterSet p)

Definition at line 20 of file MaterialBudgetForward.cc.

References book(), boundaries, constituents, detLevels, detNames, detTypes, etaRegions, edm::ParameterSet::getParameter(), testEve_cfg::level, mergeVDriftHistosByStation::name, regionTypes, and stackOrder.

                                                                     {
  
  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("MaterialBudgetForward");
  detTypes     = m_p.getParameter<std::vector<std::string> >("DetectorTypes");
  constituents = m_p.getParameter<std::vector<int> >("Constituents");
  stackOrder   = m_p.getParameter<std::vector<int> >("StackOrder");
  detNames     = m_p.getParameter<std::vector<std::string> >("DetectorNames");
  detLevels    = m_p.getParameter<std::vector<int> >("DetectorLevels");
  etaRegions   = m_p.getParameter<std::vector<double> >("EtaBoundaries");
  regionTypes  = m_p.getParameter<std::vector<int> >("RegionTypes");
  boundaries   = m_p.getParameter<std::vector<double> >("Boundaries");
  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward initialized for "
                                 << detTypes.size() << " detector types";
  unsigned int nc = 0;
  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
    edm::LogInfo("MaterialBudget") << "Type[" << ii << "] : " << detTypes[ii]
                                   << " with " << constituents[ii] <<", order "
                                   << stackOrder[ii] << " constituents --> ";
    for (int kk=0; kk<constituents[ii]; ++kk) {
      std::string name = "Unknown"; int level = -1;
      if (nc < detNames.size()) {
        name = detNames[nc]; level = detLevels[nc]; ++nc;
      }
      edm::LogInfo("MaterialBudget") << "    Constituent[" << kk << "]: "
                                     << name << " at level " << level;
    }
  }
  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward Stop condition for "
                                 << etaRegions.size() << " eta regions";
  for (unsigned int ii=0; ii<etaRegions.size(); ++ii) {
    edm::LogInfo("MaterialBudget") << "Region[" << ii << "] : Eta < " 
                                   << etaRegions[ii] << " boundary type "
                                   << regionTypes[ii] << " limit "
                                   << boundaries[ii];
  }
  book(m_p);
}
MaterialBudgetForward::~MaterialBudgetForward ( ) [virtual]

Definition at line 58 of file MaterialBudgetForward.cc.

                                              {
}
MaterialBudgetForward::MaterialBudgetForward ( const MaterialBudgetForward ) [private]

Member Function Documentation

void MaterialBudgetForward::book ( const edm::ParameterSet m_p) [private]

Definition at line 207 of file MaterialBudgetForward.cc.

References detTypes, edm::ParameterSet::getUntrackedParameter(), i, edm::Service< T >::isAvailable(), maxEta, maxSet, me100, me200, me300, me400, me500, me600, me700, me800, min, benchmark_cfg::minEta, mergeVDriftHistosByStation::name, pi, and indexGen::title.

Referenced by MaterialBudgetForward().

                                                           {

  // Book histograms
  edm::Service<TFileService> tfile;

  if ( !tfile.isAvailable() )
    throw cms::Exception("BadConfig") << "TFileService unavailable: "
                                      << "please add it to config file";

  int    binEta  = m_p.getUntrackedParameter<int>("NBinEta", 320);
  int    binPhi  = m_p.getUntrackedParameter<int>("NBinPhi", 180);
  double minEta  = m_p.getUntrackedParameter<double>("MinEta",-8.0);
  double maxEta  = m_p.getUntrackedParameter<double>("MaxEta", 8.0);
  double maxPhi  = CLHEP::pi;
  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward: Booking user "
                                 << "histos === with " << binEta << " bins "
                                 << "in eta from " << minEta << " to "
                                 << maxEta << " and " << binPhi << " bins "
                                 << "in phi from " << -maxPhi << " to "
                                 << maxPhi;

  char  name[20], title[80];
  std::string named;
  for (int i=0; i<std::min((int)(detTypes.size()+1),maxSet); i++) {
    if (i >= (int)(detTypes.size())) named = "Unknown";
    else                             named = detTypes[i];
    sprintf(name, "%d", i+100);
    sprintf(title, "MB(X0) prof Eta in %s", named.c_str());
    me100[i] =  tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
    sprintf(name, "%d", i+200);
    sprintf(title, "MB(L0) prof Eta in %s", named.c_str());
    me200[i] = tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
    sprintf(name, "%d", i+300);
    sprintf(title, "MB(Step) prof Eta in %s", named.c_str());
    me300[i] = tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
    sprintf(name, "%d", i+400);
    sprintf(title, "Eta in %s", named.c_str());
    me400[i] = tfile->make<TH1F>(name, title, binEta, minEta, maxEta);
    sprintf(name, "%d", i+500);
    sprintf(title, "MB(X0) prof Eta Phi in %s", named.c_str());
    me500[i] = tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
                                       binPhi/2, -maxPhi, maxPhi);
    sprintf(name, "%d", i+600);
    sprintf(title, "MB(L0) prof Eta Phi in %s", named.c_str());
    me600[i]= tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
                                      binPhi/2, -maxPhi, maxPhi);
    sprintf(name, "%d", i+700);
    sprintf(title, "MB(Step) prof Eta Phi in %s", named.c_str());
    me700[i]= tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
                                      binPhi/2, -maxPhi, maxPhi);
    sprintf(name, "%d", i+800);
    sprintf(title, "Eta vs Phi in %s", named.c_str());
    me800[i]= tfile->make<TH2F>(name, title, binEta/2, minEta, maxEta,
                                binPhi/2, -maxPhi, maxPhi);
  }

  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward: Booking user "
                                 << "histos done ===";

}
const MaterialBudgetForward& MaterialBudgetForward::operator= ( const MaterialBudgetForward ) [private]
bool MaterialBudgetForward::stopAfter ( const G4Step *  aStep) [private]

Definition at line 268 of file MaterialBudgetForward.cc.

References abs, boundaries, eta, etaRegions, LogDebug, regionTypes, and findQualityFiles::rr.

Referenced by update().

                                                         {

  G4ThreeVector hitPoint    = aStep->GetPreStepPoint()->GetPosition();
  if (aStep->GetPostStepPoint() != 0) 
    hitPoint = aStep->GetPostStepPoint()->GetPosition();
  double rr    = hitPoint.perp();
  double zz    = std::abs(hitPoint.z());

  bool   flag(false);
  for (unsigned int ii=0; ii<etaRegions.size(); ++ii) {
    //  LogDebug("MaterialBudget") << " MaterialBudgetForward::Eta " << eta << " in Region[" << ii << "] with " << etaRegions[ii] << " type "   << regionTypes[ii] << "|" << boundaries[ii];
    if (fabs(eta) < etaRegions[ii]) {
      if (regionTypes[ii] == 0) {
        if (rr >= boundaries[ii]-0.001) flag = true;
      } else {
        if (zz >= boundaries[ii]-0.001) flag = true;
      }
      if (flag)
        LogDebug("MaterialBudget") <<" MaterialBudgetForward::Stop after R = " 
                                   << rr << " and Z = " << zz;
      break;
    }
  }
  //  LogDebug("MaterialBudget") <<" MaterialBudgetForward:: R = " << rr << " and Z = "  << zz << " Flag " << flag;
  return flag;
}
void MaterialBudgetForward::update ( const EndOfTrack ) [private, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfTrack * >.

Definition at line 164 of file MaterialBudgetForward.cc.

References detTypes, eta, intLen, findQualityFiles::jj, LogDebug, me100, me200, me300, me400, me500, me600, me700, me800, mergeVDriftHistosByStation::name, phi, radLen, stackOrder, and stepLen.

                                                        {

  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
    for (unsigned int jj=0; jj<=detTypes.size(); ++jj) {
      if (stackOrder[jj] == (int)(ii+1)) {
        for (unsigned int kk=0; kk<=detTypes.size(); ++kk) {
          if (stackOrder[kk] == (int)(ii)) {
            radLen[jj]  += radLen[kk];
            intLen[jj]  += intLen[kk];
            LogDebug("MaterialBudget") << "MaterialBudgetForward::Add " << kk
                                       << ":" << stackOrder[kk] << " to "
                                       << jj <<":" << stackOrder[jj] <<" RadL "
                                       << radLen[kk] << " : " << radLen[jj]
                                       << " IntL " << intLen[kk] << " : "
                                       << intLen[jj] <<" StepL " << stepLen[kk]
                                       << " : " << stepLen[jj];
            break;
          }
        }
        break;
      }
    }
  }

  for (unsigned int ii=0; ii<=detTypes.size(); ++ii) {
    me100[ii]->Fill(eta, radLen[ii]);
    me200[ii]->Fill(eta, intLen[ii]);
    me300[ii]->Fill(eta, stepLen[ii]);
    me400[ii]->Fill(eta);
    me500[ii]->Fill(eta, phi, radLen[ii]);
    me600[ii]->Fill(eta, phi, intLen[ii]);
    me700[ii]->Fill(eta, phi, stepLen[ii]);
    me800[ii]->Fill(eta, phi);
    
    std::string name("Unknown");
    if (ii < detTypes.size()) name = detTypes[ii];
    LogDebug("MaterialBudget") << "MaterialBudgetForward::Volume[" << ii
                               << "]: " << name << " == Step " << stepLen[ii] 
                               << " RadL " << radLen[ii] << " IntL " 
                               << intLen[ii];
  }
}
void MaterialBudgetForward::update ( const G4Step *  ) [private, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step * >.

Definition at line 122 of file MaterialBudgetForward.cc.

References constituents, detLevels, detTypes, intLen, findQualityFiles::jj, LogDebug, logVolumes, radLen, launcher::step, stepLen, stepT, and stopAfter().

                                                      {

  //---------- each step
  G4Material * material = aStep->GetPreStepPoint()->GetMaterial();
  double step    = aStep->GetStepLength();
  double radl    = material->GetRadlen();
  double intl    = material->GetNuclearInterLength();

  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
  unsigned int indx = detTypes.size();
  unsigned int nc = 0;
  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
    for (int kk=0; kk<constituents[ii]; ++kk) {
      if (detLevels[nc+kk] <= (int)((touch->GetHistoryDepth())+1)) {
        int jj = (int)((touch->GetHistoryDepth())+1)-detLevels[nc+kk];
        if ((touch->GetVolume(jj)->GetLogicalVolume()) == logVolumes[nc+kk]) {
          indx = ii;
          break;
        }
      }
    }
    nc += (unsigned int)(constituents[ii]);
    if (indx == ii) break;
  }

  stepT         += step;
  stepLen[indx]  = stepT;
  radLen[indx]  += (step/radl);
  intLen[indx]  += (step/intl);
  LogDebug("MaterialBudget") << "MaterialBudgetForward::Step in "
                             << touch->GetVolume(0)->GetLogicalVolume()->GetName()
                             << " Index " << indx <<" Step " << step <<" RadL "
                             << step/radl << " IntL " << step/intl;

  //----- Stop tracking after selected position
  if (stopAfter(aStep)) {
    G4Track* track = aStep->GetTrack();
    track->SetTrackStatus( fStopAndKill );
  }
}
void MaterialBudgetForward::update ( const BeginOfTrack ) [private, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack * >.

Definition at line 97 of file MaterialBudgetForward.cc.

References detTypes, eta, intLen, LogDebug, phi, radLen, stepLen, and stepT.

                                                          {

  for (unsigned int ii=0; ii<(detTypes.size()+1); ++ii) {
    stepLen[ii] = 0; radLen[ii] = 0; intLen[ii] = 0;
  }

  const G4Track * aTrack = (*trk)(); // recover G4 pointer if wanted
  const G4ThreeVector& mom = aTrack->GetMomentum() ;
  if (mom.theta() != 0 ) {
    eta = mom.eta();
  } else {
    eta = -99;
  }
  phi = mom.phi();
  stepT = 0;

  double theEnergy = aTrack->GetTotalEnergy();
  int    theID     = (int)(aTrack->GetDefinition()->GetPDGEncoding());
  LogDebug("MaterialBudget") << "MaterialBudgetHcalHistos: Track "
                             << aTrack->GetTrackID() << " Code " << theID
                             << " Energy " <<theEnergy/CLHEP::GeV<<" GeV; Eta "
                             << eta << " Phi " << phi/CLHEP::deg << " PT "
                             << mom.perp()/CLHEP::GeV << " GeV *****";
}
void MaterialBudgetForward::update ( const BeginOfRun ) [private, virtual]

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfRun * >.

Definition at line 61 of file MaterialBudgetForward.cc.

References detNames, detTypes, intLen, logVolumes, mergeVDriftHistosByStation::name, radLen, stackOrder, and stepLen.

                                                     {

  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
  std::vector<G4LogicalVolume *>::const_iterator lvcite;

  unsigned int kount=detNames.size();
  for (unsigned int ii=0; ii<kount; ++ii) 
    logVolumes.push_back(0);

  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
    for (unsigned int ii=0; ii<detNames.size(); ++ii) {
      if (strcmp(detNames[ii].c_str(),(*lvcite)->GetName().c_str()) == 0) {
        logVolumes[ii] = (*lvcite);
        kount--;
        break;
      }
    }
    if (kount <= 0) break;
  }
  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward::Finds " 
                                 << detNames.size()-kount << " out of "
                                 << detNames.size() << " LV addresses";
  for (unsigned int ii=0; ii<detNames.size(); ++ii) {
    std::string name("Unknown");
    if (logVolumes[ii] != 0)  name = logVolumes[ii]->GetName();
    edm::LogInfo("MaterialBudget") << "LV[" << ii << "] : " << detNames[ii]
                                   << " Address " << logVolumes[ii] << " | " 
                                   << name;
  }

  for (unsigned int ii=0; ii<(detTypes.size()+1); ++ii) {
    stepLen.push_back(0); radLen.push_back(0); intLen.push_back(0);
  }
  stackOrder.push_back(0);
}

Member Data Documentation

std::vector<double> MaterialBudgetForward::boundaries [private]

Definition at line 50 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and stopAfter().

std::vector<int> MaterialBudgetForward::constituents [private]

Definition at line 49 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and update().

std::vector<int> MaterialBudgetForward::detLevels [private]

Definition at line 49 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and update().

std::vector<std::string> MaterialBudgetForward::detNames [private]

Definition at line 48 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and update().

std::vector<std::string> MaterialBudgetForward::detTypes [private]

Definition at line 48 of file MaterialBudgetForward.h.

Referenced by book(), MaterialBudgetForward(), and update().

double MaterialBudgetForward::eta [private]

Definition at line 58 of file MaterialBudgetForward.h.

Referenced by stopAfter(), and update().

std::vector<double> MaterialBudgetForward::etaRegions [private]

Definition at line 50 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and stopAfter().

std::vector<double> MaterialBudgetForward::intLen [private]

Definition at line 57 of file MaterialBudgetForward.h.

Referenced by update().

std::vector<G4LogicalVolume*> MaterialBudgetForward::logVolumes [private]

Definition at line 51 of file MaterialBudgetForward.h.

Referenced by update().

const int MaterialBudgetForward::maxSet = 25 [static, private]

Definition at line 52 of file MaterialBudgetForward.h.

Referenced by book().

TProfile* MaterialBudgetForward::me100[maxSet] [private]

Definition at line 55 of file MaterialBudgetForward.h.

Referenced by book(), and update().

TProfile * MaterialBudgetForward::me200[maxSet] [private]

Definition at line 55 of file MaterialBudgetForward.h.

Referenced by book(), and update().

TProfile * MaterialBudgetForward::me300[maxSet] [private]

Definition at line 55 of file MaterialBudgetForward.h.

Referenced by book(), and update().

Definition at line 53 of file MaterialBudgetForward.h.

Referenced by book(), and update().

TProfile2D* MaterialBudgetForward::me500[maxSet] [private]

Definition at line 56 of file MaterialBudgetForward.h.

Referenced by book(), and update().

TProfile2D * MaterialBudgetForward::me600[maxSet] [private]

Definition at line 56 of file MaterialBudgetForward.h.

Referenced by book(), and update().

TProfile2D * MaterialBudgetForward::me700[maxSet] [private]

Definition at line 56 of file MaterialBudgetForward.h.

Referenced by book(), and update().

Definition at line 54 of file MaterialBudgetForward.h.

Referenced by book(), and update().

double MaterialBudgetForward::phi [private]

Definition at line 58 of file MaterialBudgetForward.h.

Referenced by update().

std::vector<double> MaterialBudgetForward::radLen [private]

Definition at line 57 of file MaterialBudgetForward.h.

Referenced by update().

std::vector<int> MaterialBudgetForward::regionTypes [private]

Definition at line 49 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and stopAfter().

std::vector<int> MaterialBudgetForward::stackOrder [private]

Definition at line 49 of file MaterialBudgetForward.h.

Referenced by MaterialBudgetForward(), and update().

std::vector<double> MaterialBudgetForward::stepLen [private]

Definition at line 57 of file MaterialBudgetForward.h.

Referenced by update().

double MaterialBudgetForward::stepT [private]

Definition at line 58 of file MaterialBudgetForward.h.

Referenced by update().