ATLAS Offline Software
Loading...
Searching...
No Matches
LArWheelSolid.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "G4VGraphicsScene.hh"
6#include "G4VisExtent.hh"
7
9#include "LArWheelSolid.h"
10
11class G4NURBS;
12class G4VoxelLimits;
13class G4AffineTransform;
14
15EInside LArWheelSolid::Inside(const G4ThreeVector &inputP) const
16{
17 LWSDBG(10, std::cout << std::setprecision(25));
18 LWSDBG(1, std::cout << TypeStr() << " Inside " << MSG_VECTOR(inputP) << std::endl);
19 const EInside inside_BS = m_BoundingShape->Inside(inputP);
20 if(inside_BS == kOutside){
21 LWSDBG(2, std::cout << "outside BS" << std::endl);
22 return kOutside;
23 }
24 G4ThreeVector p( inputP );
25 int p_fan = 0;
26 const G4double d = fabs(GetCalculator()->DistanceToTheNearestFan(p, p_fan));
27 if(d > m_FHTplusT){
28 LWSDBG(2, std::cout << "outside fan d=" << d << ", m_FHTplusT=" << m_FHTplusT << std::endl);
29 return kOutside;
30 }
31 if(d < m_FHTminusT){
32 LWSDBG(2, std::cout << "inside fan d=" << d << ", m_FHTminusT=" << m_FHTminusT << ", inside_BS=" << inside(inside_BS) << std::endl);
33 return inside_BS;
34 }
35 LWSDBG(2, std::cout << "surface" << std::endl);
36 return kSurface;
37}
38
39G4ThreeVector LArWheelSolid::SurfaceNormal(const G4ThreeVector &inputP) const
40{
41 LWSDBG(1, std::cout << TypeStr() << " SurfaceNormal" << MSG_VECTOR(inputP) << std::endl);
42 EInside inside_BS = m_BoundingShape->Inside(inputP);
43 if(inside_BS != kInside){
44 LWSDBG(2, std::cout << "not inside BS" << std::endl);
45 return m_BoundingShape->SurfaceNormal(inputP);
46 }
47 G4ThreeVector p( inputP );
48 int p_fan = 0;
50 G4ThreeVector d = GetCalculator()->NearestPointOnNeutralFibre(p, p_fan);
51 d.rotateZ(inputP.phi() - p.phi()); // rotate back to initial position
52 LWSDBG(4, std::cout << "npnf" << MSG_VECTOR(d) << std::endl);
53 p = inputP - d;
54 LWSDBG(2, std::cout << "sn " << MSG_VECTOR(p.unit()) << std::endl);
55 return(p.unit());
56}
57
58G4bool LArWheelSolid::CalculateExtent(const EAxis a, const G4VoxelLimits &vl,
59 const G4AffineTransform &t, G4double &p,
60 G4double &q) const
61{
62 return m_BoundingShape->CalculateExtent(a, vl, t, p, q);
63}
64
65G4GeometryType LArWheelSolid::GetEntityType() const
66{
67 switch(m_Type){
69 return G4String("LArInnerAbsorberWheel");
70 break;
72 return G4String("LArOuterAbsorberWheel");
73 break;
75 return G4String("LArInnerElecrodWheel");
76 break;
78 return G4String("LArOuterElecrodWheel");
79 break;
81 return G4String("LArInnerAbsorberModule");
82 break;
84 return G4String("LArOuterAbsorberModule");
85 break;
87 return G4String("LArInnerElecrodModule");
88 break;
90 return G4String("LArOuterElecrodModule");
91 break;
92 case InnerGlueWheel:
93 return G4String("LArInnerGlueWheel");
94 break;
95 case OuterGlueWheel:
96 return G4String("LArOuterGlueWheel");
97 break;
98 case InnerLeadWheel:
99 return G4String("LArInnerLeadWheel");
100 break;
101 case OuterLeadWheel:
102 return G4String("LArOuterLeadWheel");
103 break;
105 return G4String("LArInnerAbsorberCone");
106 break;
108 return G4String("LArInnerElectrodCone");
109 break;
110 case InnerGlueCone:
111 return G4String("LArInnerGlueCone");
112 break;
113 case InnerLeadCone:
114 return G4String("LArInnerLeadCone");
115 break;
117 return G4String("LArOuterAbsorberFrontCone");
118 break;
120 return G4String("LArOuterElectrodFrontCone");
121 break;
123 return G4String("LArOuterGlueFrontCone");
124 break;
126 return G4String("LArOuterLeadFrontCone");
127 break;
129 return G4String("LArOuterAbsorberBackCone");
130 break;
132 return G4String("LArOuterElectrodBackCone");
133 break;
135 return G4String("LArOuterGlueBackCone");
136 break;
138 return G4String("LArOuterLeadBackCone");
139 break;
140 default:
141 G4Exception("LArWheelSolid", "UnknownSolidType", FatalException,"GetEntityType: Unknown LArWheelType.");
142 }
143 return G4String("");
144}
145
146void LArWheelSolid::DescribeYourselfTo(G4VGraphicsScene &scene) const
147{
148 scene.AddSolid(*this);
149}
150
151G4VisExtent LArWheelSolid::GetExtent() const
152{
153 return m_BoundingShape->GetExtent();
154}
155
157{
158 return m_BoundingShape->CreatePolyhedron();
159}
160
161/*
162 * returns the number of lower z boundary of z-section containing Z
163 */
164G4int LArWheelSolid::select_section(const G4double &Z) const
165{
166 for(G4int i = m_Zsect_start_search; i > 0; -- i){
167 if(Z > m_Zsect[i]) return i;
168 }
169 return 0;
170}
static Double_t a
#define LWSDBG(a, b)
@ OuterGlueWheel
@ OuterElectrodFrontCone
@ OuterLeadBackCone
@ OuterElectrodBackCone
@ OuterAbsorberFrontCone
@ OuterElectrodModule
@ OuterElectrodWheel
@ InnerAbsorberCone
@ OuterAbsorberWheel
@ OuterLeadFrontCone
@ InnerGlueWheel
@ InnerGlueCone
@ OuterGlueBackCone
@ InnerAbsorberWheel
@ InnerElectrodCone
@ InnerLeadCone
@ InnerAbsorberModule
@ OuterAbsorberModule
@ InnerLeadWheel
@ OuterLeadWheel
@ InnerElectrodModule
@ OuterAbsorberBackCone
@ OuterGlueFrontCone
@ InnerElectrodWheel
CLHEP::Hep3Vector NearestPointOnNeutralFibre(const CLHEP::Hep3Vector &p, int fan_number) const
double DistanceToTheNearestFan(CLHEP::Hep3Vector &p, int &out_fan_number) const
Determines the nearest to the input point fan.
const LArWheelSolid_t m_Type
G4int select_section(const G4double &Z) const
std::vector< G4double > m_Zsect
G4bool CalculateExtent(const EAxis, const G4VoxelLimits &, const G4AffineTransform &, G4double &, G4double &) const
G4VSolid * m_BoundingShape
void DescribeYourselfTo(G4VGraphicsScene &) const
G4int m_Zsect_start_search
G4GeometryType GetEntityType() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &) const
EInside Inside(const G4ThreeVector &) const
G4VisExtent GetExtent() const
G4Polyhedron * CreatePolyhedron() const
G4double m_FHTminusT
const LArWheelCalculator * GetCalculator(void) const
G4double m_FHTplusT