CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimG4Core/GFlash/src/GflashHadronWrapperProcess.cc

Go to the documentation of this file.
00001 //
00002 // S.Y. Jun, August 2007
00003 //
00004 #include "SimG4Core/GFlash/interface/GflashHadronWrapperProcess.h"
00005 
00006 #include "G4Track.hh"
00007 #include "G4VParticleChange.hh"
00008 #include "G4ProcessManager.hh"
00009 #include "G4ProcessVector.hh"
00010 #include "G4VProcess.hh"
00011 #include "G4GPILSelection.hh"
00012 
00013 GflashHadronWrapperProcess::GflashHadronWrapperProcess(G4String processName) :
00014   particleChange(0), 
00015   pmanager(0), 
00016   fProcessVector(0),
00017   fProcess(0) 
00018 {
00019   theProcessName = processName;
00020 }
00021 
00022 GflashHadronWrapperProcess::~GflashHadronWrapperProcess() {
00023 }
00024 
00025 G4VParticleChange* GflashHadronWrapperProcess::PostStepDoIt(const G4Track& track, const G4Step& step)
00026 {
00027   // process PostStepDoIt for the original process
00028 
00029   particleChange = pRegProcess->PostStepDoIt(track, step);
00030 
00031   // specific actions of the wrapper process 
00032 
00033   // update step/track information after PostStep of the original process is done
00034   // these steps will be repeated again without additional conflicts even if the 
00035   // parameterized physics process doesn't take over the original process
00036 
00037   particleChange->UpdateStepForPostStep(const_cast<G4Step *> (&step));
00038 
00039   // we may not want to update G4Track during this wrapper process if
00040   // G4Track accessors are not used in the G4VFastSimulationModel model
00041   //  (const_cast<G4Step *> (&step))->UpdateTrack();
00042 
00043   // update safety after each invocation of PostStepDoIts
00044   // G4double safety = std::max(endpointSafety - (endpointSafOrigin - fPostStepPoint->GetPosition()).mag(),0.);
00045   // step.GetPostStepPoint()->SetSafety(safety);
00046 
00047   // the secondaries from the original process are also created at the wrapper 
00048   // process, but they must be deleted whether the parameterized physics is  
00049   // an 'ExclusivelyForced' process or not.  Normally the secondaries will be 
00050   // created after this wrapper process in G4SteppingManager::InvokePSDIP
00051 
00052   // store the secondaries from ParticleChange to SecondaryList
00053         
00054   G4TrackVector* fSecondary = (const_cast<G4Step *> (&step))->GetfSecondary();
00055   G4int nSecondarySave = fSecondary->size();
00056 
00057   G4int num2ndaries = particleChange->GetNumberOfSecondaries();
00058   G4Track* tempSecondaryTrack;
00059 
00060   for(G4int DSecLoop=0 ; DSecLoop< num2ndaries; DSecLoop++){
00061     tempSecondaryTrack = particleChange->GetSecondary(DSecLoop);
00062     
00063     // Set the process pointer which created this track
00064     tempSecondaryTrack->SetCreatorProcess( pRegProcess );
00065 
00066     //add secondaries from this wrapper process to the existing list
00067     fSecondary->push_back( tempSecondaryTrack );
00068 
00069   } //end of loop on secondary
00070 
00071   // Now we can still impose conditions on the secondaries from ModelTrigger,
00072   // such as the number of secondaries produced by the original process as well as
00073   // those by other continuous processes - for the hadronic inelastic interaction, 
00074   // it is always true that particleChange->GetNumberOfSecondaries() > 0
00075 
00076   // at this stage, updated step information after all processes involved 
00077   // should be available
00078 
00079   //  Print(step);
00080 
00081   // call ModelTrigger for the second time - the primary loop of PostStepGPIL
00082   // is inside G4SteppingManager::DefinePhysicalStepLength().    
00083 
00084   pmanager = track.GetDefinition()->GetProcessManager();
00085 
00086   G4double testGPIL = DBL_MAX;
00087   G4double fStepLength = 0.0;
00088   G4ForceCondition fForceCondition = InActivated;
00089 
00090   fStepLength = step.GetStepLength();
00091 
00092   fProcessVector = pmanager->GetPostStepProcessVector(typeDoIt);
00093 
00094   // keep the current status of track, use fPostponeToNextEvent word for
00095   // this particular PostStep GPIL and then restore G4TrackStatus if the
00096   // paramterized physics doesn't meet trigger conditions in ModelTrigger
00097 
00098   const G4TrackStatus keepStatus = track.GetTrackStatus();
00099 
00100   (const_cast<G4Track *> (&track))->SetTrackStatus(fPostponeToNextEvent);
00101   
00102   for(G4int ipm = 0 ; ipm < fProcessVector->entries() ; ipm++) {
00103     fProcess = (*fProcessVector)(ipm);  
00104 
00105     if ( fProcess->GetProcessType() == fParameterisation ) {
00106       // test ModelTrigger via PostStepGPIL
00107 
00108       testGPIL = fProcess->PostStepGPIL(track,fStepLength,&fForceCondition );
00109 
00110       // if G4FastSimulationModel:: ModelTrigger is true, then the parameterized 
00111       // physics process takes over the current process 
00112       
00113       if( fForceCondition == ExclusivelyForced) {     
00114 
00115         // clean up memory for changing the process - counter clean up for
00116         // the secondaries created by new G4Track in G4HadronicProcess::FillTotalResult 
00117         G4int nsec = particleChange->GetNumberOfSecondaries();
00118         for(G4int DSecLoop=0 ; DSecLoop< nsec ; DSecLoop++){
00119           G4Track* tempSecondaryTrack = particleChange->GetSecondary(DSecLoop);
00120           delete tempSecondaryTrack;
00121         } 
00122         particleChange->Clear();
00123 
00124         // updating G4Step between PostStepGPIL and PostStepDoIt for the parameterized 
00125         // process may not be necessary, but do it anyway
00126 
00127         (const_cast<G4Step *> (&step))->SetStepLength(testGPIL);
00128         (const_cast<G4Track *> (&track))->SetStepLength(testGPIL);
00129 
00130         step.GetPostStepPoint()->SetStepStatus(fExclusivelyForcedProc);;
00131         step.GetPostStepPoint()->SetProcessDefinedStep(fProcess);
00132         step.GetPostStepPoint()->SetSafety(0.0);
00133 
00134         // invoke PostStepDoIt: equivalent steps for G4SteppingManager::InvokePSDIP 
00135         particleChange = fProcess->PostStepDoIt(track,step);
00136 
00137         // update PostStepPoint of Step according to ParticleChange
00138         particleChange->UpdateStepForPostStep(const_cast<G4Step *> (&step));
00139 
00140         // update G4Track according to ParticleChange after each PostStepDoIt
00141         (const_cast<G4Step *> (&step))->UpdateTrack();
00142 
00143         // update safety after each invocation of PostStepDoIts - acutally this
00144         // is not necessary for the parameterized physics process, but do it anyway
00145         step.GetPostStepPoint()->SetSafety(0.0);
00146         
00147         // additional nullification 
00148         (const_cast<G4Track *> (&track))->SetTrackStatus( particleChange->GetTrackStatus() );
00149       }
00150       else {
00151         //restore TrackStatus if fForceCondition !=  ExclusivelyForced
00152         (const_cast<G4Track *> (&track))->SetTrackStatus(keepStatus);
00153       }
00154       // assume that there is one and only one parameterized physics
00155       break;
00156     }
00157   }
00158 
00159   // remove secondaries of this wrapper process that were added to the secondary list
00160   // since they will be added in the normal stepping procedure after this->PostStepDoIt 
00161   // in G4SteppingManager::InvokePSDIP
00162 
00163   //move the iterator to the (nSecondarySave+1)th element in the secondary list
00164   G4TrackVector::iterator itv = fSecondary->begin();
00165   itv += nSecondarySave;
00166 
00167   //delete next num2ndaries tracks from the secondary list
00168   fSecondary->erase(itv,itv+num2ndaries);
00169 
00170   //end of specific actions of this wrapper process
00171 
00172   return particleChange;
00173 }
00174 
00175 void GflashHadronWrapperProcess::Print(const G4Step& step) {
00176 
00177   std::cout << " GflashHadronWrapperProcess ProcessName, PreStepPosition, preStepPoint KE, PostStepPoint KE, DeltaEnergy Nsec \n " 
00178          << step.GetPostStepPoint()->GetProcessDefinedStep()->GetProcessName() << " " 
00179          << step.GetPostStepPoint()->GetPosition() << " "  
00180          << step.GetPreStepPoint()->GetKineticEnergy()/GeV  << " "  
00181          << step.GetPostStepPoint()->GetKineticEnergy()/GeV  << " " 
00182          << step.GetDeltaEnergy()/GeV << " "
00183          << particleChange->GetNumberOfSecondaries() << std::endl;
00184 }