source: trunk/source/digits_hits/scorer/src/G4PSSphereSurfaceFlux.cc @ 1340

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

update ti head

File size: 9.2 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: G4PSSphereSurfaceFlux.cc,v 1.7 2010/07/23 04:35:38 taso Exp $
28// GEANT4 tag $Name:  $
29//
30// G4PSSphereSurfaceFlux
31#include "G4PSSphereSurfaceFlux.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 only Surface Flux.
42//  Flux version assumes only for G4Sphere shape.
43//
44// Surface is defined  at the inside of sphere.
45// Direction                  -Rmin   +Rmax
46//   0  IN || OUT            ->|<-     |
47//   1  IN                   ->|       |
48//   2  OUT                    |<-     |
49//
50// Created: 2005-11-14  Tsukasa ASO, Akinori Kimura.
51// 29-Mar-2007  T.Aso,  Bug fix for momentum direction at outgoing flux.
52// 2010-07-22   Introduce Unit specification.
53// 2010-07-22   Add weighted and divideByAre options
54//
55///////////////////////////////////////////////////////////////////////////////
56
57G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux(G4String name, 
58                                         G4int direction, G4int depth)
59    :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction),
60     weighted(true),divideByArea(true)
61{
62    DefineUnitAndCategory();
63    SetUnit("percm2");
64}
65
66G4PSSphereSurfaceFlux::G4PSSphereSurfaceFlux(G4String name, 
67                                             G4int direction,
68                                             const G4String& unit,
69                                             G4int depth)
70  :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction)
71{
72    DefineUnitAndCategory();
73    SetUnit(unit);
74}
75
76G4PSSphereSurfaceFlux::~G4PSSphereSurfaceFlux()
77{;}
78
79G4bool G4PSSphereSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
80{
81  G4StepPoint* preStep = aStep->GetPreStepPoint();
82  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
83  G4VPVParameterisation* physParam = physVol->GetParameterisation();
84  G4VSolid * solid = 0;
85  if(physParam)
86  { // for parameterized volume
87    G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
88                ->GetReplicaNumber(indexDepth);
89    solid = physParam->ComputeSolid(idx, physVol);
90    solid->ComputeDimensions(physParam,idx,physVol);
91  }
92  else
93  { // for ordinary volume
94    solid = physVol->GetLogicalVolume()->GetSolid();
95  }
96
97  G4Sphere* sphereSolid = (G4Sphere*)(solid);
98
99  G4int dirFlag =IsSelectedSurface(aStep,sphereSolid);
100  if ( dirFlag > 0 ) {
101    if ( fDirection == fFlux_InOut || fDirection == dirFlag ){
102
103      G4StepPoint* thisStep=0;
104      if ( dirFlag == fFlux_In ){
105        thisStep = preStep;
106      }else if ( dirFlag == fFlux_Out ){
107        thisStep = aStep->GetPreStepPoint();
108      }else{
109        return FALSE;
110      }
111
112      G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
113      G4ThreeVector pdirection = thisStep->GetMomentumDirection();
114      G4ThreeVector localdir  = 
115        theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
116      G4double localdirL2 = localdir.x()*localdir.x()
117        +localdir.y()*localdir.y()
118        +localdir.z()*localdir.z();
119      G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
120      G4ThreeVector localpos1 = 
121        theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
122      G4double localR2 = localpos1.x()*localpos1.x()
123        +localpos1.y()*localpos1.y()
124        +localpos1.z()*localpos1.z();
125      G4double anglefactor = (localdir.x()*localpos1.x()
126                              +localdir.y()*localpos1.y()
127                              +localdir.z()*localpos1.z())
128        /std::sqrt(localdirL2)/std::sqrt(localR2);
129
130      G4double radi   = sphereSolid->GetInsideRadius();
131      G4double dph    = sphereSolid->GetDeltaPhiAngle()/radian;
132      G4double stth   = sphereSolid->GetStartThetaAngle()/radian;
133      G4double enth   = stth+sphereSolid->GetDeltaThetaAngle()/radian;
134      G4double square = radi*radi*dph*( -std::cos(enth) + std::cos(stth) );
135
136      G4double current = 1.0;
137      if ( weighted ) thisStep->GetWeight(); // Flux (Particle Weight)
138      if ( divideByArea ) current = current/square;  // Flux with angle.
139
140      current /= anglefactor;
141
142      G4int index = GetIndex(aStep);
143      EvtMap->add(index,current);
144    }
145  }
146
147  return TRUE;
148}
149
150G4int G4PSSphereSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Sphere* sphereSolid){
151
152  G4TouchableHandle theTouchable = 
153    aStep->GetPreStepPoint()->GetTouchableHandle();
154  G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
155 
156  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
157    // Entering Geometry
158    G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
159    G4ThreeVector localpos1 = 
160      theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
161    G4double localR2 = localpos1.x()*localpos1.x()
162                      +localpos1.y()*localpos1.y()
163                      +localpos1.z()*localpos1.z();
164    //G4double InsideRadius2 =
165    //  sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
166    //if(std::fabs( localR2 - InsideRadius2 ) < kCarTolerance ){
167    G4double InsideRadius = sphereSolid->GetInsideRadius();
168    if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
169         &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
170      return fFlux_In;
171    }
172  }
173
174  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
175    // Exiting Geometry
176    G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
177    G4ThreeVector localpos2 = 
178      theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
179    G4double localR2 = localpos2.x()*localpos2.x()
180                      +localpos2.y()*localpos2.y()
181                      +localpos2.z()*localpos2.z();
182    //G4double InsideRadius2 =
183    //  sphereSolid->GetInsideRadius()*sphereSolid->GetInsideRadius();
184    //if(std::facb(localR2 - InsideRadius2) ) < kCarTolerance ){
185    G4double InsideRadius = sphereSolid->GetInsideRadius();
186    if ( localR2 > (InsideRadius-kCarTolerance)*(InsideRadius-kCarTolerance)
187         &&localR2 < (InsideRadius+kCarTolerance)*(InsideRadius+kCarTolerance)){
188      return fFlux_Out;
189    }
190  }
191
192  return -1;
193}
194
195void G4PSSphereSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
196{
197  EvtMap = new G4THitsMap<G4double>(detector->GetName(), GetName());
198  if ( HCID < 0 ) HCID = GetCollectionID(0);
199  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
200}
201
202void G4PSSphereSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
203{;}
204
205void G4PSSphereSurfaceFlux::clear(){
206  EvtMap->clear();
207}
208
209void G4PSSphereSurfaceFlux::DrawAll()
210{;}
211
212void G4PSSphereSurfaceFlux::PrintAll()
213{
214  G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
215  G4cout << " PrimitiveScorer " << GetName() <<G4endl; 
216  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
217  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
218  for(; itr != EvtMap->GetMap()->end(); itr++) {
219    G4cout << "  copy no.: " << itr->first
220           << "  current  : " << *(itr->second)/GetUnitValue()
221           << " ["<<GetUnit()<<"]"
222           << G4endl;
223  }
224}
225
226void G4PSSphereSurfaceFlux::SetUnit(const G4String& unit)
227{
228    if ( divideByArea ) {
229        CheckAndSetUnit(unit,"Per Unit Surface");
230    } else {
231        if (unit == "" ){
232            unitName = unit;
233            unitValue = 1.0;
234        }else{
235            G4String msg = "Invalid unit ["+unit+"] (Current  unit is [" +GetUnit()+"] )";
236            G4Exception(GetName(),"DetScorer0000",JustWarning,msg);
237        }
238    }
239}
240
241void G4PSSphereSurfaceFlux::DefineUnitAndCategory(){
242   // Per Unit Surface
243   new G4UnitDefinition("percentimeter2","percm2","Per Unit Surface",(1./cm2));
244   new G4UnitDefinition("permillimeter2","permm2","Per Unit Surface",(1./mm2));
245   new G4UnitDefinition("permeter2","perm2","Per Unit Surface",(1./m2));
246}
247
Note: See TracBrowser for help on using the repository browser.