Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:19:02

0001 // Tell emacs that this is a C++ source
0002 //  -*- C++ -*-.
0003 
0004 /*!
0005  * \file PHG4SectorConstructor.h
0006  * \brief Generalized sector detectors
0007  * \author Jin Huang <jhuang@bnl.gov>
0008  * \version $Revision: 1.3 $
0009  * \date $Date: 2014/05/01 19:03:26 $
0010  */
0011 
0012 #ifndef G4DETECTORS_PHG4SECTORCONSTRUCTOR_H
0013 #define G4DETECTORS_PHG4SECTORCONSTRUCTOR_H
0014 
0015 #include <Geant4/G4PhysicalConstants.hh>
0016 #include <Geant4/G4String.hh>  // for G4String
0017 #include <Geant4/G4SystemOfUnits.hh>
0018 #include <Geant4/G4Types.hh>  // for G4int
0019 
0020 class G4LogicalVolume;
0021 class G4PVPlacement;
0022 class G4VSolid;
0023 class PHG4SectorDisplayAction;
0024 class PHG4Subsystem;
0025 
0026 #include <map>
0027 #include <utility>
0028 
0029 #include <cassert>
0030 #include <cmath>
0031 #include <string>
0032 #include <vector>
0033 
0034 namespace PHG4Sector
0035 {
0036   //! layer description data
0037   //! use GEANT units!
0038   class Layer
0039   {
0040    public:
0041     Layer(  //! name base for this layer
0042         const std::string &_name,
0043 
0044         //! material name in G4
0045         const std::string &_material,
0046 
0047         //! depth in G4 units
0048         double _depth,
0049 
0050         //! percentage filled
0051         double _percentage_filled,
0052 
0053         //! whether it is active
0054         bool _active)
0055       : name(_name)
0056       , material(_material)
0057       , depth(_depth)
0058       , percentage_filled(
0059             _percentage_filled)
0060       , active(_active)
0061     {
0062     }
0063 
0064    public:
0065     //! name base for this layer
0066     std::string name;
0067 
0068     //! material name in G4
0069     std::string material;
0070 
0071     //! depth in G4 units
0072     double depth;
0073 
0074     //! percentage filled
0075     double percentage_filled;
0076 
0077     //! whether it is active
0078     bool active;
0079   };
0080 
0081   //! geometry data
0082   //! use GEANT units! Use
0083   class Sector_Geometry
0084   {
0085    public:
0086     Sector_Geometry()
0087     {
0088       SetDefault();
0089     }
0090 
0091     void
0092     SetDefault();
0093 
0094     int get_N_Sector() const
0095     {
0096       return N_Sector;
0097     }
0098 
0099     std::vector<Layer> &
0100     get_layer_list()
0101     {
0102       return layer_list;
0103     }
0104 
0105     const std::string
0106     &get_material() const
0107     {
0108       return material;
0109     }
0110 
0111     double
0112     get_max_polar_angle() const
0113     {
0114       return max_polar_angle;
0115     }
0116 
0117     double
0118     get_min_polar_angle() const
0119     {
0120       return min_polar_angle;
0121     }
0122 
0123     double
0124     get_normal_polar_angle() const
0125     {
0126       return normal_polar_angle;
0127     }
0128 
0129     double
0130     get_normal_start() const
0131     {
0132       return normal_start;
0133     }
0134 
0135     void
0136     set_N_Sector(int sector)
0137     {
0138       assert(sector >= 1);
0139 
0140       N_Sector = sector;
0141     }
0142 
0143     void
0144     set_layer_list(const std::vector<Layer> &layerList)
0145     {
0146       layer_list = layerList;
0147     }
0148 
0149     void
0150     set_material(const std::string &_material)
0151     {
0152       material = _material;
0153     }
0154 
0155     void
0156     set_max_polar_angle(double maxPolarAngle)
0157     {
0158       assert(maxPolarAngle >= 0);
0159       assert(maxPolarAngle <= pi);
0160 
0161       max_polar_angle = maxPolarAngle;
0162     }
0163 
0164     void
0165     set_min_polar_angle(double minPolarAngle)
0166     {
0167       assert(minPolarAngle >= 0);
0168       assert(minPolarAngle <= pi);
0169 
0170       min_polar_angle = minPolarAngle;
0171     }
0172 
0173     void
0174     set_normal_polar_angle(double normalPolarAngle)
0175     {
0176       assert(normalPolarAngle >= 0);
0177       assert(normalPolarAngle <= pi);
0178 
0179       normal_polar_angle = normalPolarAngle;
0180     }
0181 
0182     void
0183     set_normal_start(double normalZStart)
0184     {
0185       normal_start = normalZStart;
0186     }
0187 
0188     // derivative constants
0189 
0190     //   ! max radius from IP
0191     double
0192     get_max_R() const;
0193 
0194     double
0195     get_total_thickness() const;
0196 
0197     //! Unit
0198     static double
0199     Unit_cm()
0200     {
0201       return cm;
0202     }
0203 
0204     //! Intercept certain z point at certain polar angle
0205     void
0206     set_normal_start(const double z_intercept, const double angle_intercept)
0207     {
0208       set_normal_start(
0209           z_intercept / cos(angle_intercept) * cos(normal_polar_angle - angle_intercept));
0210     }
0211 
0212     //! Pseudorapidity
0213     static double
0214     eta_to_polar_angle(const double eta)
0215     {
0216       return 2. * atan(exp(-eta));
0217     }
0218 
0219     // layer descriptions
0220    public:
0221     typedef std::vector<Layer> t_layer_list;
0222     t_layer_list layer_list;
0223 
0224     int GetNumActiveLayers() const
0225     {
0226       int n = 0;
0227       for (t_layer_list::const_iterator it = layer_list.begin();
0228            it != layer_list.end(); ++it)
0229         if ((*it).active)
0230           n++;
0231       return n;
0232     }
0233 
0234     void
0235     AddLayer(                            //
0236         const std::string &_name,        //! name base for this layer
0237         const std::string &_material,    //! material name in G4
0238         double _depth,                   //! depth in G4 units
0239         bool _active = false,            //! active detector element for sensitive detector?
0240         double _percentage_filled = 100  //! percentage filled//
0241     )
0242     {
0243       layer_list.push_back(
0244           Layer(_name, _material, _depth, _percentage_filled, _active));
0245     }
0246 
0247     //! add Entrace window and drift volume
0248     //! Ref: P. Abbon et al. The COMPASS experiment at CERN. Nucl. Instrum. Meth., A577:455
0249     //! 518, 2007. arXiv:hep-ex/0703049, doi:10.1016/j.nima.2007.03.026. 3
0250     void
0251     AddLayers_DriftVol_COMPASS(const double drift_vol_thickness = 3 * mm);
0252 
0253     //! add HBD GEM to the layer list,
0254     //! Ref: W. Anderson et al. Design, Construction, Operation and Performance of a Hadron
0255     //! Blind Detector for the PHENIX Experiment. Nucl. Instrum. Meth., A646:35 58, 2011.
0256     //! arXiv:1103.4277, doi:10.1016/j.nima.2011.04.015. 3.5.1
0257     void
0258     AddLayers_HBD_GEM(const int n_GEM_layers = 3);
0259 
0260     //! add HBD readout,
0261     //! Ref: W. Anderson et al. Design, Construction, Operation and Performance of a Hadron
0262     //! Blind Detector for the PHENIX Experiment. Nucl. Instrum. Meth., A646:35 58, 2011.
0263     //! arXiv:1103.4277, doi:10.1016/j.nima.2011.04.015. 3.5.1
0264     void
0265     AddLayers_HBD_Readout();
0266 
0267     //! Rough AeroGel detector
0268     //! Ref: T. Iijima et al. A Novel type of proximity focusing RICH counter with multiple
0269     //! refractive index aerogel radiator. Nucl. Instrum. Meth., A548:383 390, 2005. arXiv:
0270     //! physics/0504220, doi:10.1016/j.nima.2005.05.030
0271     void
0272     AddLayers_AeroGel_ePHENIX(const double radiator_length = 2 * cm,
0273                               const double expansion_length = 18 * cm, std::string radiator = "Default");
0274 
0275    public:
0276     typedef enum
0277     {
0278 
0279       //! cone cut for the polar edge
0280       kConeEdge = 0,
0281 
0282       //! flat line edge in the azimuthal direction and along the normal_polar_angle
0283       kFlatEdge = 1
0284 
0285     } e_edge_typ;
0286 
0287     static e_edge_typ
0288     ConeEdge()
0289     {
0290       return kConeEdge;
0291     }
0292 
0293     static e_edge_typ
0294     FlatEdge()
0295     {
0296       return kFlatEdge;
0297     }
0298 
0299     e_edge_typ
0300     get_max_polar_edge() const
0301     {
0302       return max_polar_edge;
0303     }
0304 
0305     e_edge_typ
0306     get_min_polar_edge() const
0307     {
0308       return min_polar_edge;
0309     }
0310 
0311     void
0312     set_max_polar_edge(e_edge_typ maxPolarEdge)
0313     {
0314       max_polar_edge = maxPolarEdge;
0315     }
0316 
0317     void
0318     set_min_polar_edge(e_edge_typ minPolarEdge)
0319     {
0320       min_polar_edge = minPolarEdge;
0321     }
0322 
0323    private:
0324     //! number of sectors
0325     int N_Sector;
0326 
0327     //! polar angle for the normal vector
0328     double normal_polar_angle;
0329 
0330     //! polar angle for edges
0331     double min_polar_angle;
0332 
0333     //! edge type
0334     e_edge_typ min_polar_edge;
0335 
0336     //! polar angle for edges
0337     double max_polar_angle;
0338 
0339     //! edge type
0340     e_edge_typ max_polar_edge;
0341 
0342     //! distance that detector starts from the normal direction
0343     double normal_start;
0344 
0345     //! base material, usually the gas. Will fill between layers
0346     std::string material;
0347   };
0348 
0349   //! \brief Generalized detector which use sectors of flat panels to cover full azimuthal acceptance
0350   class PHG4SectorConstructor
0351   {
0352    public:
0353     PHG4SectorConstructor(const std::string &name, PHG4Subsystem *subsys);
0354     virtual ~PHG4SectorConstructor() {}
0355 
0356     void
0357     Construct_Sectors(G4LogicalVolume *WorldLog);
0358 
0359     void
0360     OverlapCheck(bool check)
0361     {
0362       overlapcheck_sector = check;
0363     }
0364 
0365     void Verbosity(int v) { m_Verbosity = v; }
0366     int Verbosity() const { return m_Verbosity; }
0367 
0368    protected:
0369     bool overlapcheck_sector;
0370 
0371     G4VSolid *
0372     Construct_Sectors_Plane(           //
0373         const std::string &name,       //
0374         const double start_z,          //
0375         const double thickness,        //
0376         G4VSolid *SecConeBoundary_Det  //
0377     );
0378 
0379    public:
0380     // properties
0381 
0382     std::string name_base;
0383 
0384     Sector_Geometry geom;
0385 
0386    private:
0387     PHG4SectorDisplayAction *m_DisplayAction;
0388     int m_Verbosity;
0389 
0390    protected:
0391     G4LogicalVolume *
0392     RegisterLogicalVolume(G4LogicalVolume *);
0393 
0394     typedef std::map<G4String, G4LogicalVolume *> map_log_vol_t;
0395     map_log_vol_t map_log_vol;
0396 
0397     G4PVPlacement *
0398     RegisterPhysicalVolume(G4PVPlacement *v, const bool active = false);
0399 
0400     typedef std::pair<G4String, G4int> phy_vol_idx_t;
0401     typedef std::map<phy_vol_idx_t, G4PVPlacement *> map_phy_vol_t;
0402     map_phy_vol_t map_phy_vol;         //! all physics volume
0403     map_phy_vol_t map_active_phy_vol;  //! active physics volume
0404   };
0405 
0406 }  // namespace PHG4Sector
0407 #endif /* PHG4SectorConstructor_H_ */