G4PSPassageCellFlux.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 // GEANT4 tag $Name: geant4-09-04 $
00029 //
00030 // G4PSPassageCellFlux
00031 #include "G4PSPassageCellFlux.hh"
00032 
00033 #include "G4SystemOfUnits.hh"
00034 #include "G4StepStatus.hh"
00035 #include "G4Track.hh"
00036 #include "G4VSolid.hh"
00037 #include "G4VPhysicalVolume.hh"
00038 #include "G4VPVParameterisation.hh"
00039 #include "G4UnitsTable.hh"
00041 // (Description)
00042 //   This is a primitive scorer class for scoring cell flux.
00043 //   The Cell Flux is defined by  a track length divided by a geometry
00044 //   volume, where only tracks passing through the geometry are taken 
00045 //  into account. e.g. the unit of Cell Flux is mm/mm3.
00046 //
00047 //   If you want to score all tracks in the geometry volume,
00048 //  please use G4PSCellFlux.
00049 //
00050 // Created: 2005-11-14  Tsukasa ASO, Akinori Kimura.
00051 // 2010-07-22   Introduce Unit specification.
00052 // 2010-07-22   Add weighted option
00053 // 
00055 
00056 G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, G4int depth)
00057   : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
00058     weighted(true)
00059 {
00060     DefineUnitAndCategory();
00061     SetUnit("percm2");
00062 }
00063 
00064 G4PSPassageCellFlux::G4PSPassageCellFlux(G4String name, const G4String& unit,
00065                                          G4int depth)
00066   : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
00067     weighted(true)
00068 {
00069     DefineUnitAndCategory();
00070     SetUnit(unit);
00071 }
00072 
00073 G4PSPassageCellFlux::~G4PSPassageCellFlux()
00074 {;}
00075 
00076 G4bool G4PSPassageCellFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
00077 {
00078 
00079   if ( IsPassed(aStep) ) {
00080     G4int idx = ((G4TouchableHistory*)
00081                  (aStep->GetPreStepPoint()->GetTouchable()))
00082                  ->GetReplicaNumber(indexDepth);
00083     G4double cubicVolume = ComputeVolume(aStep, idx);
00084 
00085     fCellFlux /= cubicVolume;
00086     G4int index = GetIndex(aStep);
00087     EvtMap->add(index,fCellFlux);
00088   }
00089 
00090   return TRUE;
00091 }
00092 
00093 G4bool G4PSPassageCellFlux::IsPassed(G4Step* aStep){
00094   G4bool Passed = FALSE;
00095 
00096   G4bool IsEnter = aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary;
00097   G4bool IsExit  = aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary;
00098 
00099   G4int  trkid  = aStep->GetTrack()->GetTrackID();
00100   G4double trklength  = aStep->GetStepLength();
00101   if ( weighted ) trklength *= aStep->GetPreStepPoint()->GetWeight();
00102 
00103   if ( IsEnter &&IsExit ){         // Passed at one step
00104     fCellFlux = trklength;         // Track length is absolutely given.
00105     Passed = TRUE;                 
00106   }else if ( IsEnter ){            // Enter a new geometry
00107     fCurrentTrkID = trkid;         // Resetting the current track.
00108     fCellFlux  = trklength;     
00109   }else if ( IsExit ){             // Exit a current geometry
00110     if ( fCurrentTrkID == trkid ) {// if the track is same as entered,
00111       fCellFlux  += trklength;     // add the track length to current one.
00112       Passed = TRUE;               
00113     }
00114   }else{                           // Inside geometry
00115     if ( fCurrentTrkID == trkid ){ // if the track is same as entered,
00116       fCellFlux  += trklength;     // adding the track length to current one.
00117     }
00118   }
00119 
00120   return Passed;
00121 }
00122 
00123 void G4PSPassageCellFlux::Initialize(G4HCofThisEvent* HCE)
00124 {
00125   fCurrentTrkID = -1;
00126 
00127   EvtMap = new G4THitsMap<G4double>(detector->GetName(),
00128                                     GetName());
00129   if ( HCID < 0 ) HCID = GetCollectionID(0);
00130   HCE->AddHitsCollection(HCID,EvtMap);
00131 
00132 }
00133 
00134 void G4PSPassageCellFlux::EndOfEvent(G4HCofThisEvent*)
00135 {;}
00136 
00137 void G4PSPassageCellFlux::clear(){
00138   EvtMap->clear();
00139 }
00140 
00141 void G4PSPassageCellFlux::DrawAll()
00142 {;}
00143 
00144 void G4PSPassageCellFlux::PrintAll()
00145 {
00146   G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
00147   G4cout << " PrimitiveScorer " << GetName() <<G4endl; 
00148   G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00149   std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00150   for(; itr != EvtMap->GetMap()->end(); itr++) {
00151     G4cout << "  copy no.: " << itr->first
00152            << "  cell flux : " << *(itr->second)/GetUnitValue()
00153            << " [" << GetUnit()
00154            << G4endl;
00155   }
00156 }
00157 
00158 void G4PSPassageCellFlux::SetUnit(const G4String& unit)
00159 {
00160     CheckAndSetUnit(unit,"Per Unit Surface");
00161 }
00162 
00163 void G4PSPassageCellFlux::DefineUnitAndCategory(){
00164    // Per Unit Surface
00165    new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
00166    new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
00167    new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
00168 }
00169 
00170 
00171 G4double G4PSPassageCellFlux::ComputeVolume(G4Step* aStep, G4int idx){
00172 
00173   G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
00174   G4VPVParameterisation* physParam = physVol->GetParameterisation();
00175   G4VSolid* solid = 0;
00176   if(physParam)
00177   { // for parameterized volume
00178     if(idx<0)
00179     {
00180       G4ExceptionDescription ED;
00181       ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
00182       G4Exception("G4PSPassageCellFlux::ComputeVolume","DetPS0013",JustWarning,ED);
00183     }
00184     solid = physParam->ComputeSolid(idx, physVol);
00185     solid->ComputeDimensions(physParam,idx,physVol);
00186   }
00187   else
00188   { // for ordinary volume
00189     solid = physVol->GetLogicalVolume()->GetSolid();
00190   }
00191   
00192   return solid->GetCubicVolume();
00193 }

Generated on Mon May 27 17:49:29 2013 for Geant4 by  doxygen 1.4.7