G4PSPassageCellFlux Class Reference

#include <G4PSPassageCellFlux.hh>

Inheritance diagram for G4PSPassageCellFlux:

G4VPrimitiveScorer G4PSPassageCellFlux3D G4PSPassageCellFluxForCylinder3D

Public Member Functions

 G4PSPassageCellFlux (G4String name, G4int depth=0)
 G4PSPassageCellFlux (G4String name, const G4String &unit, G4int depth=0)
virtual ~G4PSPassageCellFlux ()
void Weighted (G4bool flg=true)
virtual void Initialize (G4HCofThisEvent *)
virtual void EndOfEvent (G4HCofThisEvent *)
virtual void clear ()
virtual void DrawAll ()
virtual void PrintAll ()
virtual void SetUnit (const G4String &unit)

Protected Member Functions

virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *)
virtual G4bool IsPassed (G4Step *)
virtual G4double ComputeVolume (G4Step *, G4int idx)
virtual void DefineUnitAndCategory ()

Detailed Description

Definition at line 53 of file G4PSPassageCellFlux.hh.


Constructor & Destructor Documentation

G4PSPassageCellFlux::G4PSPassageCellFlux ( G4String  name,
G4int  depth = 0 
)

Definition at line 56 of file G4PSPassageCellFlux.cc.

References DefineUnitAndCategory(), and SetUnit().

00057   : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
00058     weighted(true)
00059 {
00060     DefineUnitAndCategory();
00061     SetUnit("percm2");
00062 }

G4PSPassageCellFlux::G4PSPassageCellFlux ( G4String  name,
const G4String unit,
G4int  depth = 0 
)

Definition at line 64 of file G4PSPassageCellFlux.cc.

References DefineUnitAndCategory(), and SetUnit().

00066   : G4VPrimitiveScorer(name,depth),HCID(-1),fCurrentTrkID(-1),fCellFlux(0),
00067     weighted(true)
00068 {
00069     DefineUnitAndCategory();
00070     SetUnit(unit);
00071 }

G4PSPassageCellFlux::~G4PSPassageCellFlux (  )  [virtual]

Definition at line 73 of file G4PSPassageCellFlux.cc.

00074 {;}


Member Function Documentation

void G4PSPassageCellFlux::clear (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 137 of file G4PSPassageCellFlux.cc.

00137                                {
00138   EvtMap->clear();
00139 }

G4double G4PSPassageCellFlux::ComputeVolume ( G4Step ,
G4int  idx 
) [protected, virtual]

Reimplemented in G4PSPassageCellFluxForCylinder3D.

Definition at line 171 of file G4PSPassageCellFlux.cc.

References G4VSolid::ComputeDimensions(), G4VPVParameterisation::ComputeSolid(), G4endl, G4Exception(), G4VSolid::GetCubicVolume(), G4VPhysicalVolume::GetLogicalVolume(), G4VPhysicalVolume::GetParameterisation(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPreStepPoint(), G4LogicalVolume::GetSolid(), and JustWarning.

Referenced by ProcessHits().

00171                                                                    {
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 }

void G4PSPassageCellFlux::DefineUnitAndCategory (  )  [protected, virtual]

Definition at line 163 of file G4PSPassageCellFlux.cc.

Referenced by G4PSPassageCellFlux().

00163                                                {
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 }

void G4PSPassageCellFlux::DrawAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 141 of file G4PSPassageCellFlux.cc.

00142 {;}

void G4PSPassageCellFlux::EndOfEvent ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 134 of file G4PSPassageCellFlux.cc.

00135 {;}

void G4PSPassageCellFlux::Initialize ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 123 of file G4PSPassageCellFlux.cc.

References G4HCofThisEvent::AddHitsCollection(), G4VPrimitiveScorer::detector, G4VPrimitiveScorer::GetCollectionID(), G4VPrimitiveScorer::GetName(), and G4VSensitiveDetector::GetName().

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 }

G4bool G4PSPassageCellFlux::IsPassed ( G4Step  )  [protected, virtual]

Definition at line 93 of file G4PSPassageCellFlux.cc.

References FALSE, fGeomBoundary, G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Step::GetStepLength(), G4StepPoint::GetStepStatus(), G4Step::GetTrack(), G4Track::GetTrackID(), G4StepPoint::GetWeight(), and TRUE.

Referenced by ProcessHits().

00093                                                  {
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 }

void G4PSPassageCellFlux::PrintAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 144 of file G4PSPassageCellFlux.cc.

References G4VPrimitiveScorer::detector, G4cout, G4endl, G4VPrimitiveScorer::GetName(), G4VSensitiveDetector::GetName(), G4VPrimitiveScorer::GetUnit(), and G4VPrimitiveScorer::GetUnitValue().

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 }

G4bool G4PSPassageCellFlux::ProcessHits ( G4Step ,
G4TouchableHistory  
) [protected, virtual]

Implements G4VPrimitiveScorer.

Definition at line 76 of file G4PSPassageCellFlux.cc.

References ComputeVolume(), G4VPrimitiveScorer::GetIndex(), G4Step::GetPreStepPoint(), G4StepPoint::GetTouchable(), G4VPrimitiveScorer::indexDepth, IsPassed(), and TRUE.

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 }

void G4PSPassageCellFlux::SetUnit ( const G4String unit  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 158 of file G4PSPassageCellFlux.cc.

References G4VPrimitiveScorer::CheckAndSetUnit().

Referenced by G4PSPassageCellFlux(), and G4PSPassageCellFlux3D::G4PSPassageCellFlux3D().

00159 {
00160     CheckAndSetUnit(unit,"Per Unit Surface");
00161 }

void G4PSPassageCellFlux::Weighted ( G4bool  flg = true  )  [inline]

Definition at line 62 of file G4PSPassageCellFlux.hh.

00062 { weighted = flg; }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:03 2013 for Geant4 by  doxygen 1.4.7