source: trunk/source/digits_hits/scorer/src/G4PSFlatSurfaceFlux.cc@ 1341

Last change on this file since 1341 was 1340, checked in by garnier, 15 years ago

update ti head

File size: 7.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id: G4PSFlatSurfaceFlux.cc,v 1.6 2010/07/23 04:35:38 taso Exp $
28// GEANT4 tag $Name: $
29//
30// G4PSFlatSurfaceFlux
31#include "G4PSFlatSurfaceFlux.hh"
32#include "G4StepStatus.hh"
33#include "G4Track.hh"
34#include "G4VSolid.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4VPVParameterisation.hh"
37#include "G4UnitsTable.hh"
38#include "G4GeometryTolerance.hh"
39////////////////////////////////////////////////////////////////////////////////
40// (Description)
41// This is a primitive scorer class for scoring Surface Flux.
42// Current version assumes only for G4Box shape, and the surface
43// is fixed on -Z plane.
44//
45// Surface is defined at the -Z surface.
46// Direction -Z +Z
47// 0 IN || OUT ->|<- |
48// 1 IN ->| |
49// 2 OUT |<- |
50//
51// Created: 2005-11-14 Tsukasa ASO, Akinori Kimura.
52//
53// 18-Nov-2005 T.Aso, To use always positive value for anglefactor.
54// 29-Mar-2007 T.Aso, Bug fix for momentum direction at outgoing flux.
55// 2010-07-22 Introduce Unit specification.
56// 2010-07-22 Add weighted and divideByAre options
57///////////////////////////////////////////////////////////////////////////////
58
59G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(G4String name,
60 G4int direction, G4int depth)
61 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
62 weighted(true),divideByArea(true)
63{
64 DefineUnitAndCategory();
65 SetUnit("percm2");
66}
67
68G4PSFlatSurfaceFlux::G4PSFlatSurfaceFlux(G4String name,
69 G4int direction,
70 const G4String& unit,
71 G4int depth)
72 :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction)
73{
74 DefineUnitAndCategory();
75 SetUnit(unit);
76}
77
78G4PSFlatSurfaceFlux::~G4PSFlatSurfaceFlux()
79{;}
80
81G4bool G4PSFlatSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
82{
83 G4StepPoint* preStep = aStep->GetPreStepPoint();
84 G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
85 G4VPVParameterisation* physParam = physVol->GetParameterisation();
86 G4VSolid * solid = 0;
87 if(physParam)
88 { // for parameterized volume
89 G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
90 ->GetReplicaNumber(indexDepth);
91 solid = physParam->ComputeSolid(idx, physVol);
92 solid->ComputeDimensions(physParam,idx,physVol);
93 }
94 else
95 { // for ordinary volume
96 solid = physVol->GetLogicalVolume()->GetSolid();
97 }
98
99 G4Box* boxSolid = (G4Box*)(solid);
100
101 G4int dirFlag =IsSelectedSurface(aStep,boxSolid);
102 if ( dirFlag > 0 ) {
103 if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
104
105 G4StepPoint* thisStep=0;
106 if ( dirFlag == fFlux_In ){
107 thisStep = preStep;
108 }else if ( dirFlag == fFlux_Out ){
109 thisStep = aStep->GetPostStepPoint();
110 }else{
111 return FALSE;
112 }
113
114 G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
115 G4ThreeVector pdirection = thisStep->GetMomentumDirection();
116 G4ThreeVector localdir =
117 theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
118 //
119 G4double angleFactor = localdir.z();
120 if ( angleFactor < 0 ) angleFactor *= -1.;
121 G4double flux = 1.0;
122 if ( weighted ) flux *=preStep->GetWeight(); // Current (Particle Weight)
123 //
124 G4double square = 4.*boxSolid->GetXHalfLength()*boxSolid->GetYHalfLength();
125 //
126 flux = flux/angleFactor; // Flux with angle.
127 if ( divideByArea ) flux /= square;
128 //
129 G4int index = GetIndex(aStep);
130 EvtMap->add(index,flux);
131 }
132 }
133#ifdef debug
134 G4cout << " PASSED vol "
135 << index << " trk "<<trkid<<" len " << fFlatSurfaceFlux<<G4endl;
136#endif
137
138 return TRUE;
139}
140
141G4int G4PSFlatSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Box* boxSolid){
142
143 G4TouchableHandle theTouchable =
144 aStep->GetPreStepPoint()->GetTouchableHandle();
145 G4double kCarTolerance=G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
146
147 if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
148 // Entering Geometry
149 G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
150 G4ThreeVector localpos1 =
151 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
152 if(std::fabs( localpos1.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
153 return fFlux_In;
154 }
155 }
156
157 if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
158 // Exiting Geometry
159 G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
160 G4ThreeVector localpos2 =
161 theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
162 if(std::fabs( localpos2.z() + boxSolid->GetZHalfLength())<kCarTolerance ){
163 return fFlux_Out;
164 }
165 }
166
167 return -1;
168}
169
170void G4PSFlatSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
171{
172 EvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(),
173 GetName());
174 if ( HCID < 0 ) HCID = GetCollectionID(0);
175 HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
176}
177
178void G4PSFlatSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
179{;}
180
181void G4PSFlatSurfaceFlux::clear(){
182 EvtMap->clear();
183}
184
185void G4PSFlatSurfaceFlux::DrawAll()
186{;}
187
188void G4PSFlatSurfaceFlux::PrintAll()
189{
190 G4cout << " MultiFunctionalDet " << detector->GetName() << G4endl;
191 G4cout << " PrimitiveScorer" << GetName() <<G4endl;
192 G4cout << " Number of entries " << EvtMap->entries() << G4endl;
193 std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
194 for(; itr != EvtMap->GetMap()->end(); itr++) {
195 G4cout << " copy no.: " << itr->first
196 << " flux : " << *(itr->second)/GetUnitValue()
197 << " [" << GetUnit() <<"]"
198 << G4endl;
199 }
200}
201
202void G4PSFlatSurfaceFlux::SetUnit(const G4String& unit)
203{
204 if ( divideByArea ) {
205 CheckAndSetUnit(unit,"Per Unit Surface");
206 } else {
207 if (unit == "" ){
208 unitName = unit;
209 unitValue = 1.0;
210 }else{
211 G4String msg = "Invalid unit ["+unit+"] (Current unit is [" +GetUnit()+"] )";
212 G4Exception(GetName(),"DetScorer0000",JustWarning,msg);
213 }
214 }
215}
216
217void G4PSFlatSurfaceFlux::DefineUnitAndCategory(){
218 // Per Unit Surface
219 new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
220 new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
221 new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
222}
223
224
Note: See TracBrowser for help on using the repository browser.