Class Reference for E1039 Core & Analysis Software
PHG4ConeSubsystem.cc
Go to the documentation of this file.
1 #include "PHG4ConeSubsystem.h"
2 #include "PHG4ConeDetector.h"
6 
8 #include <phool/getClass.h>
9 
10 #include <Geant4/globals.hh>
11 
12 #include <sstream>
13 
14 using namespace std;
15 
16 //_______________________________________________________________________
17 PHG4ConeSubsystem::PHG4ConeSubsystem( const std::string &name, const int lyr ):
18  PHG4Subsystem( name ),
19  detector_( 0 ),
20  steppingAction_( NULL ),
21  eventAction_(NULL),
22  place_in_x(0),
23  place_in_y(0),
24  place_in_z(0),
25  rot_in_z(0),
26  rMin1(5*cm),
27  rMax1(100*cm),
28  rMin2(5*cm),
29  rMax2(200*cm),
30  dZ(100*cm),
31  sPhi(0),
32  dPhi(2*M_PI),
33  material("Silicon"),
34  active(0),
35  layer(lyr),
36  detector_type(name),
37  superdetector("NONE")
38 {
39 
40  // put the layer into the name so we get unique names
41  // for multiple layers
42  ostringstream nam;
43  nam << name << "_" << lyr;
44  Name(nam.str().c_str());
45 }
46 
47 //_______________________________________________________________________
49 {
50  PHNodeIterator iter( topNode );
51  PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode","DST" ));
52 
53  // create detector
54  detector_ = new PHG4ConeDetector(topNode, Name(), layer);
55  detector_->SetR1(rMin1, rMax1);
56  detector_->SetR2(rMin2, rMax2);
57  detector_->SetZlength(dZ);
58  detector_->SetPhi(sPhi, dPhi);
59  detector_->SetPlace(place_in_x, place_in_y, place_in_z);
60  detector_->SetZRot(rot_in_z);
61  detector_->SetMaterial(material);
62  detector_->SetActive(active);
63  detector_->SuperDetector(superdetector);
64  detector_->OverlapCheck(overlapcheck);
65  if (active)
66  {
67  ostringstream nodename;
68  if (superdetector != "NONE")
69  {
70  nodename << "G4HIT_" << superdetector;
71  }
72  else
73  {
74  nodename << "G4HIT_" << detector_type << "_" << layer;
75  }
76  // create hit list
77  PHG4HitContainer* block_hits = findNode::getClass<PHG4HitContainer>( topNode ,nodename.str().c_str());
78  if ( !block_hits )
79  {
80 
81  dstNode->addNode( new PHIODataNode<PHObject>( block_hits = new PHG4HitContainer(nodename.str()),nodename.str().c_str(),"PHObject" ));
82 
83  }
84  // create stepping action
85  steppingAction_ = new PHG4ConeSteppingAction(detector_);
86 
87  eventAction_ = new PHG4EventActionClearZeroEdep(topNode, nodename.str());
88  }
89  return 0;
90 
91 }
92 
93 //_______________________________________________________________________
94 int
96 {
97  // pass top node to stepping action so that it gets
98  // relevant nodes needed internally
99  if (steppingAction_)
100  {
101  steppingAction_->SetInterfacePointers( topNode );
102  }
103  return 0;
104 }
105 
106 
107 //_______________________________________________________________________
109 { return detector_; }
110 
111 //_______________________________________________________________________
113 { return steppingAction_; }
114 
115 
116 //_______________________________________________________________________
117 void
118 PHG4ConeSubsystem::Set_eta_range(G4double etaMin, G4double etaMax)
119 {
120  G4double thetaMin = 2*atan(exp(-etaMax));
121  G4double thetaMax = 2*atan(exp(-etaMin));
122 
123  G4double z1 = place_in_z - dZ;
124  G4double z2 = place_in_z + dZ;
125 
126  rMin1 = z1*tan(thetaMin);
127  rMax1 = z1*tan(thetaMax);
128 
129  rMin2 = z2*tan(thetaMin);
130  rMax2 = z2*tan(thetaMax);
131 }
#define NULL
Definition: Pdb.h:9
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
PHBoolean addNode(PHNode *)
void SetZlength(const G4double a)
set length in Z
void SetActive(const int i=1)
void SetR2(const G4double min, const G4double max)
set inner and outter radius2
void SetMaterial(const std::string &mat)
void SuperDetector(const std::string &name)
void SetPhi(const G4double a, const G4double b)
set phi offset and extention
void SetZRot(const G4double z_angle)
void SetPlace(const G4double place_x, const G4double place_y, const G4double place_z)
void SetR1(const G4double min, const G4double max)
set inner and outter radius1
virtual void SetInterfacePointers(PHCompositeNode *)
reimplemented from base class
int process_event(PHCompositeNode *)
event processing
virtual PHG4Detector * GetDetector(void) const
accessors (reimplemented)
virtual PHG4SteppingAction * GetSteppingAction(void) const
return pointer to this subsystem stepping action
PHG4ConeSubsystem(const std::string &name="CONE", const int layer=0)
constructor
int Init(PHCompositeNode *)
init
void Set_eta_range(G4double etaMin, G4double etaMax)
set rmaximum and minimums according to the eta range
base class for phenix detector creation
Definition: PHG4Detector.h:14
virtual void OverlapCheck(const bool chk=true)
Definition: PHG4Detector.h:53
PHNode * findFirst(const std::string &, const std::string &)