G4PSCellFlux Class Reference

#include <G4PSCellFlux.hh>

Inheritance diagram for G4PSCellFlux:

G4VPrimitiveScorer G4PSCellFlux3D G4PSCellFluxForCylinder3D

Public Member Functions

 G4PSCellFlux (G4String name, G4int depth=0)
 G4PSCellFlux (G4String name, const G4String &unit, G4int depth=0)
virtual ~G4PSCellFlux ()
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 G4double ComputeVolume (G4Step *, G4int idx)
virtual void DefineUnitAndCategory ()

Detailed Description

Definition at line 58 of file G4PSCellFlux.hh.


Constructor & Destructor Documentation

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

Definition at line 56 of file G4PSCellFlux.cc.

References DefineUnitAndCategory(), and SetUnit().

00057     :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
00058 {
00059     DefineUnitAndCategory();
00060     SetUnit("percm2");
00061     //verboseLevel = 10;
00062 }

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

Definition at line 64 of file G4PSCellFlux.cc.

References DefineUnitAndCategory(), and SetUnit().

00065     :G4VPrimitiveScorer(name,depth),HCID(-1),weighted(true)
00066 {
00067     DefineUnitAndCategory();
00068     SetUnit(unit);
00069 }

G4PSCellFlux::~G4PSCellFlux (  )  [virtual]

Definition at line 71 of file G4PSCellFlux.cc.

00072 {;}


Member Function Documentation

void G4PSCellFlux::clear (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 103 of file G4PSCellFlux.cc.

00103                         {
00104   EvtMap->clear();
00105 }

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

Reimplemented in G4PSCellFluxForCylinder3D.

Definition at line 136 of file G4PSCellFlux.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().

00136                                                             {
00137 
00138   G4VPhysicalVolume* physVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
00139   G4VPVParameterisation* physParam = physVol->GetParameterisation();
00140   G4VSolid* solid = 0;
00141   if(physParam)
00142   { // for parameterized volume
00143     if(idx<0)
00144     {
00145       G4ExceptionDescription ED;
00146       ED << "Incorrect replica number --- GetReplicaNumber : " << idx << G4endl;
00147       G4Exception("G4PSCellFlux::ComputeVolume","DetPS0001",JustWarning,ED);
00148     }
00149     solid = physParam->ComputeSolid(idx, physVol);
00150     solid->ComputeDimensions(physParam,idx,physVol);
00151   }
00152   else
00153   { // for ordinary volume
00154     solid = physVol->GetLogicalVolume()->GetSolid();
00155   }
00156   
00157   return solid->GetCubicVolume();
00158 }

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

Definition at line 129 of file G4PSCellFlux.cc.

Referenced by G4PSCellFlux().

00129                                         {
00130    // Per Unit Surface
00131    new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
00132    new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
00133    new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
00134 }

void G4PSCellFlux::DrawAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 107 of file G4PSCellFlux.cc.

00108 {;}

void G4PSCellFlux::EndOfEvent ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 100 of file G4PSCellFlux.cc.

00101 {;}

void G4PSCellFlux::Initialize ( G4HCofThisEvent  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 92 of file G4PSCellFlux.cc.

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

00093 {
00094   EvtMap = new G4THitsMap<G4double>(detector->GetName(),
00095                                     GetName());
00096   if ( HCID < 0 ) HCID = GetCollectionID(0);
00097   HCE->AddHitsCollection(HCID,EvtMap);
00098 }

void G4PSCellFlux::PrintAll (  )  [virtual]

Reimplemented from G4VPrimitiveScorer.

Definition at line 110 of file G4PSCellFlux.cc.

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

00111 {
00112   G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
00113   G4cout << " PrimitiveScorer " << GetName() <<G4endl; 
00114   G4cout << " Number of entries " << EvtMap->entries() << G4endl;
00115   std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
00116   for(; itr != EvtMap->GetMap()->end(); itr++) {
00117     G4cout << "  copy no.: " << itr->first
00118            << "  cell flux : " << *(itr->second)/GetUnitValue() 
00119            << " [" << GetUnit() << "]"
00120            << G4endl;
00121   }
00122 }

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

Implements G4VPrimitiveScorer.

Definition at line 74 of file G4PSCellFlux.cc.

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

00075 {
00076   G4double stepLength = aStep->GetStepLength();
00077   if ( stepLength == 0. ) return FALSE;
00078 
00079   G4int idx = ((G4TouchableHistory*)
00080                (aStep->GetPreStepPoint()->GetTouchable()))
00081                ->GetReplicaNumber(indexDepth);
00082   G4double cubicVolume = ComputeVolume(aStep, idx);
00083 
00084   G4double CellFlux = stepLength / cubicVolume;
00085   if (weighted) CellFlux *= aStep->GetPreStepPoint()->GetWeight(); 
00086   G4int index = GetIndex(aStep);
00087   EvtMap->add(index,CellFlux);
00088 
00089   return TRUE;
00090 }

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

Reimplemented from G4VPrimitiveScorer.

Definition at line 124 of file G4PSCellFlux.cc.

References G4VPrimitiveScorer::CheckAndSetUnit().

Referenced by G4PSCellFlux(), and G4PSCellFlux3D::G4PSCellFlux3D().

00125 {
00126     CheckAndSetUnit(unit,"Per Unit Surface");
00127 }

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

Definition at line 65 of file G4PSCellFlux.hh.

00065 { weighted = flg; }


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