ATLAS Offline Software
Loading...
Searching...
No Matches
LArWheelSolid.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef GEO2G4_LARWHEELSOLID_H
6#define GEO2G4_LARWHEELSOLID_H
7#ifndef PORTABLE_LAR_SHAPE
10#endif
11#include "G4VSolid.hh"
12
13// set this to allow debug output in LArWheelSolid methods
14// disabled by default to avoid any performance degradation
15//#define DEBUG_LARWHEELSOLID
16
17// set this to use native G4 FanBound's methods for DisToIn
18// instead of local calculations
19//#define LARWHEELSOLID_USE_FANBOUND
20
21// set this to use BoundingShape's methods for DisToOut
22// instead of local calculations
23#define LARWHEELSOLID_USE_BS_DTO
24
25// change this to have more z sections
26#define LARWHEELSOLID_ZSECT_MULT 1
27
28
29// set this to check in dti and dto functions if particle direction
30// pointing inside or outside of volume to return zero fast when it is required by spec.
31// currently at development stage, requires accurate surface normal calculations
32//#define CHECK_DIRTONORM_ANGLE_ON_SURFACE
33
34#ifdef DEBUG_LARWHEELSOLID
35#define LWSDBG(a, b) if(Verbose >= a) b
36#define MSG_VECTOR(v) "(" << v.x() << ", " << v.y() << ", " << v.z() << ")"
37//#define LWS_HARD_TEST_DTI
38//#define LWS_HARD_TEST_DTO
39#else
40#define LWSDBG(a, b)
41#endif
42
43// Forward declarations.
44class G4VGraphicsScene;
45class G4VisExtent;
46class G4Polyhedron;
47class G4NURBS;
48class G4VoxelLimits;
49class G4AffineTransform;
50class G4Polycone;
52class TF1;
53class LArFanSections;
54class G4Polyhedra;
55struct EMECData;
56
57#include "LArWheelSolid_type.h"
58
60{
61 switch(type){
62 case InnerAbsorberWheel: return("InnerAbsorberWheel");
63 case OuterAbsorberWheel: return("OuterAbsorberWheel");
64 case InnerElectrodWheel: return("InnerElectrodWheel");
65 case OuterElectrodWheel: return("OuterElectrodWheel");
66 case InnerAbsorberModule: return("InnerAbsorberModule");
67 case OuterAbsorberModule: return("OuterAbsorberModule");
68 case InnerElectrodModule: return("InnerElectrodModule");
69 case OuterElectrodModule: return("OuterElectrodModule");
70 case InnerLeadWheel: return("InnerLeadWheel");
71 case OuterLeadWheel: return("OuterLeadWheel");
72 case InnerGlueWheel: return("InnerGlueWheel");
73 case OuterGlueWheel: return("OuterGlueWheel");
74 case InnerAbsorberCone: return("InnerAbsorberCone");
75 case InnerElectrodCone: return("InnerElectrodCone");
76 case InnerGlueCone: return("InnerGlueCone");
77 case InnerLeadCone: return("InnerLeadCone");
78 case OuterAbsorberFrontCone: return("OuterAbsorberFrontCone");
79 case OuterElectrodFrontCone: return("OuterElectrodFrontCone");
80 case OuterGlueFrontCone: return("OuterGlueFrontCone");
81 case OuterLeadFrontCone: return("OuterLeadFrontCone");
82 case OuterAbsorberBackCone: return("OuterAbsorberBackCone");
83 case OuterElectrodBackCone: return("OuterElectrodBackCone");
84 case OuterGlueBackCone: return("OuterGlueBackCone");
85 case OuterLeadBackCone: return("OuterLeadBackCone");
86 }
87 return("unknown");
88}
89
90class ATLAS_NOT_THREAD_SAFE LArWheelSolid : public G4VSolid
91#ifndef PORTABLE_LAR_SHAPE
92, public AthMessaging
93#endif
94{
95public:
96
97 LArWheelSolid(const G4String& name, LArWheelSolid_t type,
98 G4int zside = 1,
99 LArWheelCalculator *calc = 0,
100 const EMECData *emecData=0
101 );
102 virtual ~LArWheelSolid();
103
104 // Mandatory for custom solid Geant4 functions
105 EInside Inside(const G4ThreeVector&) const;
106 G4double DistanceToIn(const G4ThreeVector&,
107 const G4ThreeVector&) const;
108 G4double DistanceToIn(const G4ThreeVector&) const;
109 G4double DistanceToOut(const G4ThreeVector&,
110 const G4ThreeVector&,
111 const G4bool calcNorm = false,
112 G4bool* validNorm = 0,
113 G4ThreeVector* n = 0) const;
114 G4double DistanceToOut(const G4ThreeVector&) const;
115 G4ThreeVector SurfaceNormal (const G4ThreeVector&) const;
116 G4bool CalculateExtent(const EAxis,
117 const G4VoxelLimits&,
118 const G4AffineTransform&,
119 G4double&,
120 G4double&) const;
121 G4GeometryType GetEntityType() const;
122 void DescribeYourselfTo(G4VGraphicsScene&) const;
123 G4VisExtent GetExtent() const;
124 G4Polyhedron* CreatePolyhedron() const;
125
126 // 07-Feb-2003 WGS: For compatibility with Geant 4.5.0
127 virtual std::ostream& StreamInfo(std::ostream& os) const { return os; }
128
129 const G4VSolid *GetBoundingShape(void) const { return m_BoundingShape; }
130 const LArWheelCalculator *GetCalculator(void) const { return m_Calculator; }
131 LArWheelSolid_t GetType(void) const { return m_Type; }
132
133#ifndef PORTABLE_LAR_SHAPE
134 G4ThreeVector GetPointOnSurface(void) const;
135 G4double GetCubicVolume(void);
136 G4double GetSurfaceArea(void);
137#endif
138
139private:
140 static const G4double s_Tolerance;
141 static const G4double s_AngularTolerance;
142 static const G4double s_IterationPrecision;
143 static const G4double s_IterationPrecision2;
144 static const unsigned int s_IterationsLimit;
145
146 G4bool m_IsOuter;
151 G4double m_MinPhi;
152 G4double m_MaxPhi;
153 const G4double m_PhiPosition;
155#ifdef LARWHEELSOLID_USE_FANBOUND
156 G4VSolid* m_FanBound;
157#endif
158
159 std::vector<G4double> m_Zsect;
161
163
164 // z at outer wheel "bend"
165 G4double m_Zmid;
166 // Special limit, used in dto
167 G4double m_Ymin;
168 // limits for use in service functions
170 //artificial level to distinguish between inner and outer cones
171 G4double m_Ymid;
172
173 void inner_solid_init(const G4String &);
174 void outer_solid_init(const G4String &);
175 void set_phi_size(void);
176
177 virtual G4double distance_to_in(G4ThreeVector &,
178 const G4ThreeVector &, int) const;
179 G4double in_iteration_process(const G4ThreeVector &,
180 G4double, G4ThreeVector &, int) const;
181 G4double search_for_nearest_point(const G4ThreeVector &, const G4double,
182 const G4ThreeVector &, int) const;
183 G4bool search_for_most_remoted_point(const G4ThreeVector &,
184 const G4ThreeVector &,
185 G4ThreeVector &, const int) const;
186 G4double out_iteration_process(const G4ThreeVector &,
187 G4ThreeVector &, const int) const;
188
189 typedef enum {
192 } FanBoundExit_t;
193
194 FanBoundExit_t find_exit_point(const G4ThreeVector &p,
195 const G4ThreeVector &v,
196 G4ThreeVector &q) const;
197 G4bool fs_cross_lower(const G4ThreeVector &p, const G4ThreeVector &v,
198 G4ThreeVector &q) const;
199 G4bool fs_cross_upper(const G4ThreeVector &p, const G4ThreeVector &v,
200 G4ThreeVector &q) const;
201 G4bool check_D(G4double &b,
202 G4double A, G4double B, G4double C, G4bool) const;
203
204 G4int select_section(const G4double &Z) const;
205
206#ifndef PORTABLE_LAR_SHAPE
207 EInside Inside_accordion(const G4ThreeVector&) const;
208 void get_point_on_accordion_surface(G4ThreeVector &) const;
209 void get_point_on_polycone_surface(G4ThreeVector &) const;
210 void get_point_on_flat_surface(G4ThreeVector &) const;
211 void set_failover_point(G4ThreeVector &p, const char *m = 0) const;
212
213 G4double get_area_on_polycone(void) const;
214 G4double get_area_on_face(void) const;
215 G4double get_area_on_side(void) const;
216
217 G4double get_area_at_r(G4double r) const;
218 G4double get_length_at_r(G4double r) const;
219
220 void test(void);
221 void clean_tests(void);
222 void init_tests(void);
223#endif
224
225protected:
226
228
230 friend double LArWheelSolid_fcn_area(double *, double *);
231 friend double LArWheelSolid_fcn_vol(double *, double *);
232 friend double LArWheelSolid_fcn_area_on_pc(double *, double *);
233 friend double LArWheelSolid_get_dl(double *, double *, G4int);
234 friend double LArWheelSolid_fcn_side_area(double *, double *);
235
236#ifdef DEBUG_LARWHEELSOLID
237 static const char* inside(EInside i)
238 {
239 switch(i){
240 case kInside: return "inside"; break;
241 case kSurface: return "surface"; break;
242 case kOutside: return "outside"; break;
243 }
244 return "unknown";
245 }
246
247public:
248 static G4int Verbose;
249 void SetVerbose(G4int v){ Verbose = v; }
250#ifndef PORTABLE_LAR_SHAPE
251 G4bool test_dti(const G4ThreeVector &p,
252 const G4ThreeVector &v, const G4double distance) const;
253 G4bool test_dto(const G4ThreeVector &p,
254 const G4ThreeVector &v, const G4double distance) const;
255#endif
256private:
257 const char *TypeStr(void) const { return LArWheelSolidTypeString(m_Type); }
258#endif
259};
260
261#endif // GEO2G4_LARWHEELSOLID_H
262
double LArWheelSolid_get_dl(double *x, double *par, G4int side)
double LArWheelSolid_fcn_area(double *x, double *p)
double LArWheelSolid_fcn_vol(double *x, double *p)
double LArWheelSolid_fcn_side_area(double *x, double *p)
double LArWheelSolid_fcn_area_on_pc(double *x, double *p)
const char * LArWheelSolidTypeString(LArWheelSolid_t type)
LArWheelSolid_t
@ OuterGlueWheel
@ OuterElectrodFrontCone
@ OuterLeadBackCone
@ OuterElectrodBackCone
@ OuterAbsorberFrontCone
@ OuterElectrodModule
@ OuterElectrodWheel
@ InnerAbsorberCone
@ OuterAbsorberWheel
@ OuterLeadFrontCone
@ InnerGlueWheel
@ InnerGlueCone
@ OuterGlueBackCone
@ InnerAbsorberWheel
@ InnerElectrodCone
@ InnerLeadCone
@ InnerAbsorberModule
@ OuterAbsorberModule
@ InnerLeadWheel
@ OuterLeadWheel
@ InnerElectrodModule
@ OuterAbsorberBackCone
@ OuterGlueFrontCone
@ InnerElectrodWheel
Define macros for attributes used to control the static checker.
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
This class separates some of the geometry details of the LAr endcap.
LArWheelCalculator * m_Calculator
G4double out_iteration_process(const G4ThreeVector &, G4ThreeVector &, const int) const
void outer_solid_init(const G4String &)
static const G4double s_IterationPrecision2
const LArWheelSolid_t m_Type
static const G4double s_IterationPrecision
static const G4double s_Tolerance
LArWheelSolid_t GetType(void) const
LArWheelSolid(const G4String &name, LArWheelSolid_t type, G4int zside=1, LArWheelCalculator *calc=0, const EMECData *emecData=0)
G4double in_iteration_process(const G4ThreeVector &, G4double, G4ThreeVector &, int) const
std::vector< G4double > m_Zsect
G4bool CalculateExtent(const EAxis, const G4VoxelLimits &, const G4AffineTransform &, G4double &, G4double &) const
G4double DistanceToIn(const G4ThreeVector &, const G4ThreeVector &) const
const G4VSolid * GetBoundingShape(void) const
G4double m_MinPhi
G4VSolid * m_BoundingShape
void DescribeYourselfTo(G4VGraphicsScene &) const
static const unsigned int s_IterationsLimit
G4int m_Zsect_start_search
static const G4double s_AngularTolerance
G4GeometryType GetEntityType() const
virtual G4double distance_to_in(G4ThreeVector &, const G4ThreeVector &, int) const
G4double DistanceToOut(const G4ThreeVector &, const G4ThreeVector &, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &) const
G4double m_FanPhiAmplitude
virtual std::ostream & StreamInfo(std::ostream &os) const
G4double m_MaxPhi
EInside Inside(const G4ThreeVector &) const
G4VisExtent GetExtent() const
G4Polyhedron * CreatePolyhedron() const
const G4double m_PhiPosition
LArFanSections * m_fs
G4bool search_for_most_remoted_point(const G4ThreeVector &, const G4ThreeVector &, G4ThreeVector &, const int) const
void inner_solid_init(const G4String &)
G4double search_for_nearest_point(const G4ThreeVector &, const G4double, const G4ThreeVector &, int) const
G4double m_FHTminusT
G4double m_FanHalfThickness
const LArWheelCalculator * GetCalculator(void) const
G4double m_FHTplusT
int r
Definition globals.cxx:22
struct color C
@ Verbose
Definition ZDCMsg.h:18
hold the test vectors and ease the comparison