source: trunk/source/digits_hits/scorer/src/G4PSCylinderSurfaceFlux.cc @ 1051

Last change on this file since 1051 was 1012, checked in by garnier, 15 years ago

update

File size: 7.4 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: G4PSCylinderSurfaceFlux.cc,v 1.3 2008/12/28 21:10:51 asaim Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// // G4PSCylinderSurfaceFlux
31#include "G4PSCylinderSurfaceFlux.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 G4Tubs shape, and the surface
43//  is fixed on inner plane of the tube.
44//
45// Surface is defined at the innner surface of the tube.
46// Direction                   R    R+dR
47//   0  IN || OUT            ->|<-  |
48//   1  IN                   ->|    |
49//   2  OUT                    |<-  |
50//
51// Created: 2007-03-29  Tsukasa ASO
52///////////////////////////////////////////////////////////////////////////////
53
54G4PSCylinderSurfaceFlux::G4PSCylinderSurfaceFlux(G4String name, 
55                                                 G4int direction, G4int depth)
56  :G4VPrimitiveScorer(name,depth),HCID(-1),fDirection(direction)
57{;}
58
59G4PSCylinderSurfaceFlux::~G4PSCylinderSurfaceFlux()
60{;}
61
62G4bool G4PSCylinderSurfaceFlux::ProcessHits(G4Step* aStep,G4TouchableHistory*)
63{
64  G4StepPoint* preStep = aStep->GetPreStepPoint();
65  G4StepPoint* postStep = aStep->GetPreStepPoint();
66
67  G4VPhysicalVolume* physVol = preStep->GetPhysicalVolume();
68  G4VPVParameterisation* physParam = physVol->GetParameterisation();
69  G4VSolid * solid = 0;
70  if(physParam)
71  { // for parameterized volume
72    G4int idx = ((G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable()))
73                ->GetReplicaNumber(indexDepth);
74    solid = physParam->ComputeSolid(idx, physVol);
75    solid->ComputeDimensions(physParam,idx,physVol);
76  }
77  else
78  { // for ordinary volume
79    solid = physVol->GetLogicalVolume()->GetSolid();
80  }
81
82//  if( solid->GetEntityType() != "G4Tubs" ){
83//    G4Exception("G4PSCylinderSurfaceFluxScorer. - Solid type is not supported.");
84//    return FALSE;
85//  }
86  G4Tubs* tubsSolid = (G4Tubs*)(solid);
87 
88  G4int dirFlag =IsSelectedSurface(aStep,tubsSolid);
89 
90  if ( dirFlag > 0 ){
91    if (fDirection == fFlux_InOut || dirFlag == fDirection ){
92
93      G4StepPoint* thisStep=0;
94      if ( dirFlag == fFlux_In ){
95        thisStep = preStep;
96      }else if ( dirFlag == fFlux_Out ){
97        thisStep = postStep;
98      }else{
99        return FALSE;
100      }
101 
102      G4TouchableHandle theTouchable = thisStep->GetTouchableHandle();
103      G4ThreeVector pdirection = thisStep->GetMomentumDirection();
104      G4ThreeVector localdir  = 
105        theTouchable->GetHistory()->GetTopTransform().TransformAxis(pdirection);
106      G4ThreeVector position = thisStep->GetPosition();
107      G4ThreeVector localpos  =
108        theTouchable->GetHistory()->GetTopTransform().TransformAxis(position);
109      G4double angleFactor = (localdir.x()*localpos.x()+localdir.y()*localpos.y())
110        /std::sqrt(localdir.x()*localdir.x()
111                   +localdir.y()*localdir.y()+localdir.z()*localdir.z())
112        /std::sqrt(localpos.x()*localpos.x()+localpos.y()*localpos.y());
113   
114      if ( angleFactor < 0 ) angleFactor *= -1.;
115      G4double square = 2.*tubsSolid->GetZHalfLength()
116        *tubsSolid->GetInnerRadius()* tubsSolid->GetDeltaPhiAngle()/radian;
117   
118      G4double flux = preStep->GetWeight(); 
119      // Current (Particle Weight)
120
121      flux = flux/angleFactor/square;   
122      //Flux with angle.
123      G4int index = GetIndex(aStep);
124      EvtMap->add(index,flux);
125      return TRUE;
126    }else{
127      return FALSE;
128    }
129  }else{
130      return FALSE;
131  }
132}
133
134G4int G4PSCylinderSurfaceFlux::IsSelectedSurface(G4Step* aStep, G4Tubs* tubsSolid){
135
136  G4TouchableHandle theTouchable = 
137    aStep->GetPreStepPoint()->GetTouchableHandle();
138  G4double kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance(); 
139
140  if (aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary ){
141    // Entering Geometry
142    G4ThreeVector stppos1= aStep->GetPreStepPoint()->GetPosition();
143    G4ThreeVector localpos1 = 
144      theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos1);
145    if ( std::fabs(localpos1.z()) > tubsSolid->GetZHalfLength() ) return -1;
146    if(std::fabs( localpos1.x()*localpos1.x()+localpos1.y()*localpos1.y() ) 
147       - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius())<kCarTolerance ){
148      return fFlux_In;
149    }
150  }
151
152  if (aStep->GetPostStepPoint()->GetStepStatus() == fGeomBoundary ){
153    // Exiting Geometry
154    G4ThreeVector stppos2= aStep->GetPostStepPoint()->GetPosition();
155    G4ThreeVector localpos2 = 
156      theTouchable->GetHistory()->GetTopTransform().TransformPoint(stppos2);
157    if ( std::fabs(localpos2.z()) > tubsSolid->GetZHalfLength() ) return -1;
158    if(std::fabs( localpos2.x()*localpos2.x()+localpos2.y()*localpos2.y() ) 
159       - (tubsSolid->GetInnerRadius()*tubsSolid->GetInnerRadius())<kCarTolerance ){
160      return fFlux_Out;
161    }
162  }
163
164  return -1;
165}
166
167void G4PSCylinderSurfaceFlux::Initialize(G4HCofThisEvent* HCE)
168{
169  EvtMap = new G4THitsMap<G4double>(GetMultiFunctionalDetector()->GetName(),
170                                    GetName());
171  if ( HCID < 0 ) HCID = GetCollectionID(0);
172  HCE->AddHitsCollection(HCID, (G4VHitsCollection*)EvtMap);
173}
174
175void G4PSCylinderSurfaceFlux::EndOfEvent(G4HCofThisEvent*)
176{;}
177
178void G4PSCylinderSurfaceFlux::clear(){
179  EvtMap->clear();
180}
181
182void G4PSCylinderSurfaceFlux::DrawAll()
183{;}
184
185void G4PSCylinderSurfaceFlux::PrintAll()
186{
187  G4cout << " MultiFunctionalDet  " << detector->GetName() << G4endl;
188  G4cout << " PrimitiveScorer" << GetName() <<G4endl; 
189  G4cout << " Number of entries " << EvtMap->entries() << G4endl;
190  std::map<G4int,G4double*>::iterator itr = EvtMap->GetMap()->begin();
191  for(; itr != EvtMap->GetMap()->end(); itr++) {
192    G4cout << "  copy no.: " << itr->first
193           << "  flux  : " << *(itr->second)*cm*cm << " [cm^-2]"
194           << G4endl;
195  }
196}
197
Note: See TracBrowser for help on using the repository browser.