CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

GflashEMShowerModel Class Reference

#include <GflashEMShowerModel.h>

List of all members.

Public Member Functions

void DoIt (const G4FastTrack &, G4FastStep &)
 GflashEMShowerModel (G4String name, G4Envelope *env, edm::ParameterSet parSet)
G4bool IsApplicable (const G4ParticleDefinition &)
G4bool ModelTrigger (const G4FastTrack &)
 ~GflashEMShowerModel ()

Private Member Functions

G4bool excludeDetectorRegion (const G4FastTrack &fastTrack)
void makeHits (const G4FastTrack &fastTrack)
void updateGflashStep (G4ThreeVector position, G4double time)

Private Attributes

G4Navigator * theGflashNavigator
G4Step * theGflashStep
G4TouchableHandle theGflashTouchableHandle
edm::ParameterSet theParSet
GflashEMShowerProfiletheProfile
bool theWatcherOn

Detailed Description

Definition at line 55 of file GflashEMShowerModel.h.


Constructor & Destructor Documentation

GflashEMShowerModel::GflashEMShowerModel ( G4String  name,
G4Envelope *  env,
edm::ParameterSet  parSet 
)

Definition at line 22 of file GflashEMShowerModel.cc.

References edm::ParameterSet::getParameter(), theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, and theWatcherOn.

  : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {

  theWatcherOn = parSet.getParameter<bool>("watcherOn");

  theProfile = new GflashEMShowerProfile(parSet);

  theGflashStep = new G4Step();
  theGflashTouchableHandle = new G4TouchableHistory();
  theGflashNavigator = new G4Navigator();

}
GflashEMShowerModel::~GflashEMShowerModel ( )

Definition at line 37 of file GflashEMShowerModel.cc.

References theGflashStep, and theProfile.

                                          {

  if(theProfile) delete theProfile;
  if(theGflashStep) delete theGflashStep;
}

Member Function Documentation

void GflashEMShowerModel::DoIt ( const G4FastTrack &  fastTrack,
G4FastStep &  fastStep 
)

Definition at line 88 of file GflashEMShowerModel.cc.

References DeDxDiscriminatorTools::charge(), relval_parameters_module::energy, Gflash::findShowerType(), GflashEMShowerProfile::initialize(), makeHits(), GflashEMShowerProfile::parameterization(), position, and theProfile.

                                                                                 {

  // Kill the parameterised particle:
  fastStep.KillPrimaryTrack();
  fastStep.ProposePrimaryTrackPathLength(0.0);

  //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
  G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
  G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
  G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
  G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
  G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
  G4int showerType = Gflash::findShowerType(position);

  // Do actual parameterization. The result of parameterization is gflashHitList
  theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
  theProfile->parameterization();

  //make hits
  makeHits(fastTrack);
}
G4bool GflashEMShowerModel::excludeDetectorRegion ( const G4FastTrack &  fastTrack) [private]

Definition at line 165 of file GflashEMShowerModel.cc.

References eta().

Referenced by ModelTrigger().

                                                                              {

  G4bool isExcluded=false;

  //exclude regions where geometry are complicated
  //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5 
  G4double eta =   fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
  if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) return true;

  return isExcluded;
}
G4bool GflashEMShowerModel::IsApplicable ( const G4ParticleDefinition &  particleType)

Definition at line 43 of file GflashEMShowerModel.cc.

                                                                                 { 

  return ( &particleType == G4Electron::ElectronDefinition() ||
           &particleType == G4Positron::PositronDefinition() ); 

}
void GflashEMShowerModel::makeHits ( const G4FastTrack &  fastTrack) [private]

Definition at line 110 of file GflashEMShowerModel.cc.

References GflashEMShowerProfile::getGflashHitList(), SteppingAction::m_g4StepSignal, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, theWatcherOn, and updateGflashStep().

Referenced by DoIt().

                                                               {

  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();

  theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));

  theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
    (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
  theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());

  std::vector<GflashHit>::const_iterator spotIter    = gflashHitList.begin();
  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();

  for( ; spotIter != spotIterEnd; spotIter++){

      //put touchable for each hit so that touchable history keeps track of each step.
      theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
                                                                    theGflashTouchableHandle, false);
      updateGflashStep(spotIter->getPosition(),spotIter->getTime());

      // if there is a watcher defined in a job and the flag is turned on
      if(theWatcherOn) {
        SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
        userSteppingAction->m_g4StepSignal(theGflashStep);
      }

      // Send G4Step information to Hit/Digi if the volume is sensitive
      // Copied from G4SteppingManager.cc
    
      G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
      if( aCurrentVolume == 0 ) continue;

      G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
      if(lv->GetRegion()->GetName() != "CaloRegion") continue;

      theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
      G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
      
      if( aSensitive == 0 ) continue;

      theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
      aSensitive->Hit(theGflashStep);
  }

}
G4bool GflashEMShowerModel::ModelTrigger ( const G4FastTrack &  fastTrack)

Definition at line 51 of file GflashEMShowerModel.cc.

References excludeDetectorRegion().

                                                                       {

  // Mininum energy cutoff to parameterize
  if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < 1.0*GeV) return false;
  if(excludeDetectorRegion(fastTrack)) return false;

  G4bool trigger = fastTrack.GetPrimaryTrack()->GetDefinition() == G4Electron::ElectronDefinition() || 
    fastTrack.GetPrimaryTrack()->GetDefinition() == G4Positron::PositronDefinition();

  if(!trigger) return false;

  // This will be changed accordingly when the way dealing with CaloRegion changes later.
  G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
  G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
  if( pCurrentVolume == 0) return false;

  G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
  if(lv->GetRegion()->GetName() != "CaloRegion") return false;

  // The parameterization starts inside crystals
  std::size_t pos1 = lv->GetName().find("EBRY");
  std::size_t pos2 = lv->GetName().find("EFRY");
  /*
  std::size_t pos3 = lv->GetName().find("HVQ");
  std::size_t pos4 = lv->GetName().find("HF");
  if(pos1 == std::string::npos && pos2 == std::string::npos &&
     pos3 == std::string::npos && pos4 == std::string::npos) return false;
  */
  //@@@for now, HF is not a part of Gflash Envelopes
  if(pos1 == std::string::npos && pos2 == std::string::npos ) return false;

  return true;

}
void GflashEMShowerModel::updateGflashStep ( G4ThreeVector  position,
G4double  time 
) [private]

Definition at line 156 of file GflashEMShowerModel.cc.

References theGflashStep, and theGflashTouchableHandle.

Referenced by makeHits().

{
  theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
  theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
  theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
  theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
}

Member Data Documentation

Definition at line 78 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().

G4TouchableHandle GflashEMShowerModel::theGflashTouchableHandle [private]

Definition at line 79 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), makeHits(), and updateGflashStep().

Definition at line 72 of file GflashEMShowerModel.h.

Definition at line 75 of file GflashEMShowerModel.h.

Referenced by DoIt(), GflashEMShowerModel(), makeHits(), and ~GflashEMShowerModel().

Definition at line 73 of file GflashEMShowerModel.h.

Referenced by GflashEMShowerModel(), and makeHits().