Class Reference for E1039 Core & Analysis Software
SQG4DipoleMagnetDetector.cc
Go to the documentation of this file.
2 
3 #include <phparameter/PHParameters.h>
4 #include <g4main/PHG4Utils.h>
6 #include <phool/PHIODataNode.h>
7 #include <phool/getClass.h>
8 #include <db_svc/DbSvc.h>
9 
10 #include <TString.h>
11 #include <TSQLStatement.h>
12 
13 #include <Geant4/G4Box.hh>
14 #include <Geant4/G4Tubs.hh>
15 #include <Geant4/G4SubtractionSolid.hh>
16 #include <Geant4/G4Colour.hh>
17 #include <Geant4/G4LogicalVolume.hh>
18 #include <Geant4/G4Material.hh>
19 #include <Geant4/G4PVPlacement.hh>
20 #include <Geant4/G4SystemOfUnits.hh>
21 #include <Geant4/G4UserLimits.hh>
22 #include <Geant4/G4VisAttributes.hh>
23 #include <Geant4/G4NistManager.hh>
24 
25 #include <iostream>
26 #include <sstream>
27 #include <map>
28 
29 SQG4DipoleMagnetDetector::SQG4DipoleMagnetDetector(PHCompositeNode* node, PHParameters* parameters, const std::string& dnam, const int lyr):
30  PHG4Detector(node, dnam),
31  params(parameters),
32  block_physi(nullptr),
33  layer(lyr)
34 {}
35 
37 
38 bool SQG4DipoleMagnetDetector::IsInBlock(G4VPhysicalVolume* volume) const
39 {
40  if(volume == block_physi)
41  {
42  return true;
43  }
44  return false;
45 }
46 
47 void SQG4DipoleMagnetDetector::Construct(G4LogicalVolume* logicWorld)
48 {
49  std::string dbfile = params->get_string_param("geomdb");
50  std::string vName = params->get_string_param("magname");
51  std::cout << "SQDipoleMagnet - begin construction of " << vName << " from " << dbfile << std::endl;
52 
53  DbSvc* db_svc = new DbSvc(DbSvc::LITE, dbfile.c_str());
54  std::unique_ptr<TSQLStatement> stmt;
55 
56  // Verbosity(10);
57  // overlapcheck = true;
58 
59  //Elements and materials are constructed in PHG4Reco::DefineMaterials() - we start from shape contruction
60 
61  //Load all shapes
62  std::map<int, G4VSolid*> solids;
63 
65  stmt.reset(db_svc->Process( "SELECT sID, sName, xLength, yLength, zLength FROM SolidBoxes WHERE sName like '" + vName + "%'"));
66  while(stmt->NextResultRow() && (!stmt->IsNull(0)))
67  {
68  int sID = stmt->GetInt(0);
69  if(solids.find(sID) != solids.end())
70  {
71  std::cout << "ERROR - SQDipoleMagnet Construction: duplicate solid ID " << sID << std::endl;
72  return;
73  }
74 
75  G4String sName = stmt->GetString(1);
76  double xLength = stmt->GetDouble(2)*cm;
77  double yLength = stmt->GetDouble(3)*cm;
78  double zLength = stmt->GetDouble(4)*cm;
79  if(Verbosity() > 1)
80  {
81  std::cout << "SQDipoleMagnet Construct: create solid box " << sID << " " << sName << " " << xLength << " " << yLength << " " << zLength << std::endl;
82  }
83  solids[sID] = new G4Box(sName, xLength/2., yLength/2., zLength/2.);
84  }
85 
87  stmt.reset(db_svc->Process( "SELECT sID, sName, length, radiusMin, radiusMax FROM SolidTubes WHERE sName like '" + vName + "%'"));
88  while(stmt->NextResultRow() && (!stmt->IsNull(0)))
89  {
90  int sID = stmt->GetInt(0);
91  if(solids.find(sID) != solids.end())
92  {
93  std::cout << "ERROR - SQDipoleMagnet Construction: duplicate solid ID " << sID << std::endl;
94  return;
95  }
96 
97  G4String sName = stmt->GetString(1);
98  double length = stmt->GetDouble(2)*cm;
99  double radiusMin = stmt->GetDouble(3)*cm;
100  double radiusMax = stmt->GetDouble(4)*cm;
101  if(Verbosity() > 1)
102  {
103  std::cout << "SQDipoleMagnet Construct: create solid tube " << sID << " " << sName << " " << length << " " << radiusMax << " " << radiusMin << std::endl;
104  }
105  solids[sID] = new G4Tubs(sName, radiusMin, radiusMax, length/2., 0., 360.*deg);
106  }
107 
109  stmt.reset(db_svc->Process( "SELECT sID, sName, shellID, holeID, rotX, rotY, rotZ, posX, posY, posZ FROM SubtractionSolids WHERE sName like '" + vName + "%'"));
110  while(stmt->NextResultRow() && (!stmt->IsNull(0)))
111  {
112  int sID = stmt->GetInt(0);
113  if(solids.find(sID) != solids.end())
114  {
115  std::cout << "ERROR - SQDipoleMagnet Construction: duplicate solid ID " << sID << std::endl;
116  return;
117  }
118 
119  G4String sName = stmt->GetString(1);
120  int shellID = stmt->GetInt(2);
121  int holeID = stmt->GetInt(3);
122  if(Verbosity() > 1)
123  {
124  std::cout << "SQDipoleMagnet Construct: create solid subtraction " << sID << " " << sName << " " << shellID << " " << holeID << std::endl;
125  }
126  if(solids.find(holeID) == solids.end() || solids.find(shellID) == solids.end())
127  {
128  std::cout << "ERROR - SQDipoleMagnet Construction: cannot find solid component for " << sName << std::endl;
129  return;
130  }
131 
132  G4RotationMatrix rot;
133  rot.rotateX(stmt->GetDouble(4)*rad);
134  rot.rotateY(stmt->GetDouble(5)*rad);
135  rot.rotateZ(stmt->GetDouble(6)*rad);
136 
137  G4ThreeVector pos;
138  pos.setX(stmt->GetDouble(7)*cm);
139  pos.setY(stmt->GetDouble(8)*cm);
140  pos.setZ(stmt->GetDouble(9)*cm);
141 
142  G4Transform3D trans(rot, pos);
143  solids[sID] = new G4SubtractionSolid(sName, solids[shellID], solids[holeID], trans);
144  }
145 
146  // Create logical volumes using Table LogicalVolumes
147  std::map<int, G4LogicalVolume*> logicals;
148  stmt.reset(db_svc->Process( "SELECT lvID, lvName, sID, mName FROM LogicalVolumes WHERE lvName like '" + vName + "%'"));
149  while(stmt->NextResultRow() && (!stmt->IsNull(0)))
150  {
151  int lvID = stmt->GetInt(0);
152  if(logicals.find(lvID) != logicals.end())
153  {
154  std::cout << "ERROR - SQDipoleMagnet Construction: duplicate logical volume ID " << lvID << std::endl;
155  return;
156  }
157 
158  G4String lvName = stmt->GetString(1);
159  int sID = stmt->GetInt(2);
160  G4String mName = stmt->GetString(3);
161  if(Verbosity() > 1)
162  {
163  std::cout << "SQDipoleMagnet Construct: create logical volume " << lvID << " " << lvName << " " << sID << " " << mName << std::endl;
164  }
165  logicals[lvID] = new G4LogicalVolume(solids[sID], G4Material::GetMaterial(mName), lvName);
166  }
167 
168  // Initialize reference counter for each logical volume
169  std::map<int, int> lvRefCount;
170  for(auto it = logicals.begin(); it != logicals.end(); ++it)
171  {
172  lvRefCount[it->first] = 0;
173  }
174 
175  // Create container volume first by looking for the entry in PhysicalVolumes with motherID=-1
176  int nTopPV = 0;
177  int topLVID = -1;
178  int topPVID = -1;
179  std::string topPVName;
180  stmt.reset(db_svc->Process( "SELECT pvID, pvName, lvID FROM PhysicalVolumes WHERE motherID=0 AND depth=1 AND pvName like '" + vName + "%'"));
181  while(stmt->NextResultRow())
182  {
183  ++nTopPV;
184 
185  topPVID = stmt->GetInt(0);
186  topPVName = stmt->GetString(1);
187  topLVID = stmt->GetInt(2);
188  }
189 
190  if(nTopPV != 1)
191  {
192  std::cout << "ERROR - SQDipoleMagnet Construction: top volume is ambiguous " << nTopPV << std::endl;
193  return;
194  }
195 
196  if(logicals.find(topLVID) == logicals.end())
197  {
198  std::cout << "ERROR - SQDipoleMagnet Construction: top logical volume is missing " << std::endl;
199  return;
200  }
201 
202  if(Verbosity() > 1)
203  {
204  std::cout << "SQDipoleMagnet Construct: create top physical volume " << topPVID << " " << topLVID << " " << topPVName << std::endl;
205  }
206 
207  // Now create the top volume
208  G4RotationMatrix topRot;
209  topRot.rotateX(params->get_double_param("rot_x")*rad);
210  topRot.rotateY(params->get_double_param("rot_y")*rad);
211  topRot.rotateZ(params->get_double_param("rot_z")*rad);
212 
213  G4ThreeVector topLoc;
214  topLoc.setX(params->get_double_param("place_x")*cm);
215  topLoc.setY(params->get_double_param("place_y")*cm);
216  topLoc.setZ(params->get_double_param("place_z")*cm);
217 
218  block_physi = new G4PVPlacement(G4Transform3D(topRot, topLoc), logicals[topLVID], topPVName.c_str(), logicWorld, false, lvRefCount[topLVID], overlapcheck);
219  ++lvRefCount[topLVID];
220 
221  // The rest of the physical volumes
222  stmt.reset(db_svc->Process( "SELECT pvName, lvID, motherID, xRel, yRel, zRel, rotX, rotY, rotZ FROM PhysicalVolumes WHERE motherID>0 AND pvName like '" + vName + "%'"));
223  while(stmt->NextResultRow() && (!stmt->IsNull(0)))
224  {
225  //may need to check duplicate pvID?
226  G4String pvName = stmt->GetString(0);
227  int lvID = stmt->GetInt(1);
228  int motherID = stmt->GetInt(2);
229 
230  G4ThreeVector pos;
231  pos.setX(stmt->GetDouble(3)*cm);
232  pos.setY(stmt->GetDouble(4)*cm);
233  pos.setZ(stmt->GetDouble(5)*cm);
234 
235  G4RotationMatrix rot;
236  rot.rotateX(stmt->GetDouble(6)*rad);
237  rot.rotateY(stmt->GetDouble(7)*rad);
238  rot.rotateZ(stmt->GetDouble(8)*rad);
239 
240  if(Verbosity() > 1)
241  {
242  std::cout << "SQDipoleMagnet Construct: create physical volume " << pvName << " " << lvID << " " << motherID << std::endl;
243  }
244 
245  new G4PVPlacement(G4Transform3D(rot, pos), logicals[lvID], pvName.c_str(), logicals[motherID], false, lvRefCount[lvID], overlapcheck);
246  ++lvRefCount[lvID];
247  }
248 }
Standard interface with SQL database.
Definition: DbSvc.h:15
TSQLStatement * Process(const char *query)
Definition: DbSvc.cc:176
@ LITE
Definition: DbSvc.h:17
base class for phenix detector creation
Definition: PHG4Detector.h:14
virtual int Verbosity() const
Definition: PHG4Detector.h:47
bool overlapcheck
Definition: PHG4Detector.h:62
double get_double_param(const std::string &name) const
std::string get_string_param(const std::string &name) const
bool IsInBlock(G4VPhysicalVolume *) const
virtual void Construct(G4LogicalVolume *world)
construct
virtual ~SQG4DipoleMagnetDetector(void)
destructor
SQG4DipoleMagnetDetector(PHCompositeNode *node, PHParameters *parameters, const std::string &dnam="BLOCK", const int lyr=0)
constructor