Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4tgbVolume.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4tgbVolume.cc 69803 2013-05-15 15:24:50Z gcosmo $
28 //
29 //
30 // class G4tgbVolume
31 
32 // History:
33 // - Created. P.Arce, CIEMAT (November 2007)
34 // -------------------------------------------------------------------------
35 
36 #include "G4tgbVolume.hh"
37 
38 #include "G4PhysicalConstants.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4tgbVolumeMgr.hh"
41 #include "G4tgbMaterialMgr.hh"
43 #include "G4tgbPlaceParamLinear.hh"
44 #include "G4tgbPlaceParamSquare.hh"
45 #include "G4tgbPlaceParamCircle.hh"
46 
47 #include "G4tgrSolid.hh"
48 #include "G4tgrSolidBoolean.hh"
49 #include "G4tgrVolume.hh"
50 #include "G4tgrVolumeDivision.hh"
51 #include "G4tgrVolumeAssembly.hh"
52 #include "G4tgrVolumeMgr.hh"
53 #include "G4tgrPlace.hh"
54 #include "G4tgrPlaceSimple.hh"
55 #include "G4tgrPlaceDivRep.hh"
57 #include "G4tgrUtils.hh"
58 
59 #include "G4VSolid.hh"
60 #include "G4UnionSolid.hh"
61 #include "G4SubtractionSolid.hh"
62 #include "G4IntersectionSolid.hh"
63 #include "G4LogicalVolume.hh"
64 #include "G4VPhysicalVolume.hh"
65 #include "G4PVPlacement.hh"
66 #include "G4PVDivision.hh"
67 #include "G4PVReplica.hh"
68 #include "G4PVParameterised.hh"
69 #include "G4Box.hh"
70 #include "G4Tubs.hh"
71 #include "G4Cons.hh"
72 #include "G4Trap.hh"
73 #include "G4Sphere.hh"
74 #include "G4Orb.hh"
75 #include "G4Trd.hh"
76 #include "G4Para.hh"
77 #include "G4Torus.hh"
78 #include "G4Hype.hh"
79 #include "G4Polycone.hh"
80 #include "G4Polyhedra.hh"
81 #include "G4EllipticalTube.hh"
82 #include "G4Ellipsoid.hh"
83 #include "G4EllipticalCone.hh"
84 #include "G4Hype.hh"
85 #include "G4Tet.hh"
86 #include "G4TwistedBox.hh"
87 #include "G4TwistedTrap.hh"
88 #include "G4TwistedTrd.hh"
89 #include "G4TwistedTubs.hh"
90 #include "G4AssemblyVolume.hh"
91 #include "G4BREPSolidBox.hh"
92 #include "G4BREPSolidCylinder.hh"
93 #include "G4BREPSolidCone.hh"
94 #include "G4BREPSolidSphere.hh"
95 #include "G4BREPSolidTorus.hh"
96 #include "G4BREPSolidPCone.hh"
97 #include "G4BREPSolidPolyhedra.hh"
98 #include "G4BREPSolidOpenPCone.hh"
99 #include "G4TessellatedSolid.hh"
100 #include "G4TriangularFacet.hh"
101 #include "G4QuadrangularFacet.hh"
102 #include "G4ExtrudedSolid.hh"
103 
104 #include "G4VisExtent.hh"
105 #include "G4Material.hh"
106 #include "G4RotationMatrix.hh"
107 #include "G4ReflectionFactory.hh"
108 
109 #include "G4VisAttributes.hh"
110 #include "G4RegionStore.hh"
111 #include "G4tgrMessenger.hh"
112 #include "G4UIcommand.hh"
113 #include "G4GeometryTolerance.hh"
114 
115 //-------------------------------------------------------------------
117  : theTgrVolume(0), theG4AssemblyVolume(0)
118 {
119 }
120 
121 
122 //-------------------------------------------------------------------
124 {
125 }
126 
127 
128 //-------------------------------------------------------------------
130 {
131  theTgrVolume = vol;
132  theG4AssemblyVolume = 0;
133 }
134 
135 
136 //-------------------------------------------------------------------
138  const G4LogicalVolume* parentLV )
139 {
140 #ifdef G4VERBOSE
142  {
143  G4cout << G4endl << "@@@ G4tgbVolume::ConstructG4Volumes - " << GetName() << G4endl;
144  if( place && parentLV ) G4cout << " place in LV " << parentLV->GetName() << G4endl;
145  }
146 #endif
148  G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( GetName() );
149  G4bool bFirstCopy = false;
150  if( logvol == 0 )
151  {
152  bFirstCopy = true;
153  if( theTgrVolume->GetType() != "VOLDivision" )
154  {
155  //--- If first time build solid and LogVol
156  G4VSolid* solid = FindOrConstructG4Solid( theTgrVolume->GetSolid() );
157  if( solid != 0 ) // for G4AssemblyVolume it is 0
158  {
159  g4vmgr->RegisterMe( solid );
160  logvol = ConstructG4LogVol( solid );
161  g4vmgr->RegisterMe( logvol );
162  g4vmgr->RegisterChildParentLVs( logvol, parentLV );
163  }
164  }
165  else
166  {
167  return;
168  }
169  }
170  //--- Construct PhysVol
171  G4VPhysicalVolume* physvol = ConstructG4PhysVol( place, logvol, parentLV );
172  if( physvol != 0 ) // 0 for G4AssemblyVolumes
173  {
174  g4vmgr->RegisterMe( physvol );
175 
176  if( logvol == 0 ) // case of divisions
177  {
178  logvol = physvol->GetLogicalVolume();
179  }
180  }
181  else
182  {
183  return;
184  }
185 
186  //--- If first copy build children placements in this LogVol
187  if(bFirstCopy)
188  {
189  std::pair<G4mmapspl::iterator, G4mmapspl::iterator> children
191  G4mmapspl::iterator cite;
192  for( cite = children.first; cite != children.second; cite++ )
193  {
194  //----- Call G4tgrPlace ->constructG4Volumes
195  //---- find G4tgbVolume corresponding to the G4tgrVolume
196  // pointed by G4tgrPlace
197  G4tgrPlace* pl = const_cast<G4tgrPlace*>((*cite).second);
198  G4tgbVolume* svol = g4vmgr->FindVolume( pl->GetVolume()->GetName() );
199  //--- find copyNo
200 #ifdef G4VERBOSE
202  {
203  G4cout << " G4tgbVolume::ConstructG4Volumes - construct daughter " << pl->GetVolume()->GetName() << " # " << pl->GetCopyNo() << G4endl;
204  }
205 #endif
206  svol->ConstructG4Volumes( pl, logvol );
207  }
208  }
209 
210 }
211 
212 
213 
214 //-------------------------------------------------------------------
216 {
217  G4double angularTolerance = G4GeometryTolerance::GetInstance()
219 
220  if( sol == 0 ) { return 0; }
221 
222 #ifdef G4VERBOSE
224  {
225  G4cout << " G4tgbVolume::FindOrConstructG4Solid():" << G4endl
226  << " SOLID = " << sol << G4endl
227  << " " << sol->GetName() << " of type " << sol->GetType()
228  << G4endl;
229  }
230 #endif
231 
232  //----- Check if solid exists already
234  ->FindG4Solid( sol->GetName() );
235  if( solid ) { return solid; }
236 
237  // Give 'sol' as Boolean solids needs to call this method twice
238 
239 #ifdef G4VERBOSE
241  {
242  G4cout << " G4tgbVolume::FindOrConstructG4Solid() - "
243  << sol->GetSolidParams().size() << G4endl;
244  }
245 #endif
246 
247  std::vector<G4double> solParam;
248 
249  // In case of BOOLEAN solids, solidParams are taken from components
250 
251  if( sol->GetSolidParams().size() == 1)
252  {
253  solParam = * sol->GetSolidParams()[ 0 ];
254  }
255 
256  //----------- instantiate the appropiate G4VSolid type
257  G4String stype = sol->GetType();
258  G4String sname = sol->GetName();
259 
260  if( stype == "BOX" )
261  {
262  CheckNoSolidParams( stype, 3, solParam.size() );
263  solid = new G4Box( sname, solParam[0], solParam[1], solParam[2] );
264 
265  }
266  else if( stype == "TUBE" )
267  {
268  CheckNoSolidParams( stype, 3, solParam.size() );
269  solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2],
270  0.*deg, 360.*deg );
271  }
272  else if( stype == "TUBS" )
273  {
274  CheckNoSolidParams( stype, 5, solParam.size() );
275  G4double phiDelta = solParam[4];
276  if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
277  solid = new G4Tubs( sname, solParam[0], solParam[1], solParam[2],
278  solParam[3], phiDelta );
279  }
280  else if( stype == "TRAP" )
281  {
282  if( solParam.size() == 11 )
283  {
284  solid = new G4Trap( sname, solParam[0], solParam[1], solParam[2],
285  solParam[3], solParam[4], solParam[5], solParam[6],
286  solParam[7], solParam[8], solParam[9], solParam[10] );
287  }
288  else if( solParam.size() == 4 )
289  {
290  solid = new G4Trap( sname, solParam[0], solParam[1]/deg,
291  solParam[2]/deg, solParam[3]);
292  }
293  else
294  {
295  G4String ErrMessage1 = "Solid type " + stype;
296  G4String ErrMessage2 = " should have 11 or 4 parameters,\n";
297  G4String ErrMessage3 = "and it has "
298  + G4UIcommand::ConvertToString(G4int(solParam.size()));
299  G4String ErrMessage = ErrMessage1 + ErrMessage2 + ErrMessage3 + " !";
300  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
301  "InvalidSetup", FatalException, ErrMessage);
302  return 0;
303  }
304 
305  }
306  else if( stype == "TRD" )
307  {
308  CheckNoSolidParams( stype, 5, solParam.size() );
309  solid = new G4Trd( sname, solParam[0], solParam[1], solParam[2],
310  solParam[3], solParam[4] );
311  }
312  else if( stype == "PARA" )
313  {
314  CheckNoSolidParams( stype, 6, solParam.size() );
315  solid = new G4Para( sname, solParam[0], solParam[1], solParam[2],
316  solParam[3], solParam[4], solParam[5] );
317  }
318  else if( stype == "CONE" )
319  {
320  CheckNoSolidParams( stype, 5, solParam.size() );
321  solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2],
322  solParam[3], solParam[4], 0., 360.*deg);
323  }
324  else if( stype == "CONS" )
325  {
326  CheckNoSolidParams( stype, 7, solParam.size() );
327  G4double phiDelta = solParam[6];
328  if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
329  solid = new G4Cons( sname, solParam[0], solParam[1], solParam[2],
330  solParam[3], solParam[4], solParam[5], phiDelta);
331  }
332  else if( stype == "SPHERE" )
333  {
334  CheckNoSolidParams( stype, 6, solParam.size() );
335  G4double phiDelta = solParam[3];
336  if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
337  G4double thetaDelta = solParam[5];
338  if( std::fabs(thetaDelta - pi) < angularTolerance ) { thetaDelta = pi; }
339  solid = new G4Sphere( sname, solParam[0], solParam[1], solParam[2],
340  phiDelta, solParam[4], thetaDelta);
341  }
342  else if( stype == "ORB" )
343  {
344  CheckNoSolidParams( stype, 1, solParam.size() );
345  solid = new G4Orb( sname, solParam[0] );
346  }
347  else if( stype == "TORUS" )
348  {
349  CheckNoSolidParams( stype, 5, solParam.size() );
350  G4double phiDelta = solParam[4];
351  if( std::fabs(phiDelta - twopi) < angularTolerance ) { phiDelta = twopi; }
352  solid = new G4Torus( sname, solParam[0], solParam[1], solParam[2],
353  solParam[3], phiDelta );
354  }
355  else if( stype == "POLYCONE" )
356  {
357  size_t nplanes = size_t(solParam[2]);
358  G4bool genericPoly = false;
359  if( solParam.size() == 3+nplanes*3 )
360  {
361  genericPoly = true;
362  }
363  else if( solParam.size() == 3+nplanes*2 )
364  {
365  genericPoly = false;
366  }
367  else
368  {
369  G4String Err1 = "Solid type " + stype + " should have ";
370  G4String Err2 = G4UIcommand::ConvertToString(G4int(3+nplanes*3))
371  + " (Z,Rmin,Rmax)\n";
372  G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(3+nplanes*2));
373  G4String Err4 = " (RZ corners) parameters,\n";
374  G4String Err5 = "and it has "
375  + G4UIcommand::ConvertToString(G4int(solParam.size()));
376  G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !";
377  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
378  "InvalidSetup", FatalException, ErrMessage);
379  return 0;
380  }
381 
382  if( genericPoly )
383  {
384  std::vector<G4double>* z_p = new std::vector<G4double>;
385  std::vector<G4double>* rmin_p = new std::vector<G4double>;
386  std::vector<G4double>* rmax_p = new std::vector<G4double>;
387  for( size_t ii = 0; ii < nplanes; ii++ )
388  {
389  (*z_p).push_back( solParam[3+3*ii] );
390  (*rmin_p).push_back( solParam[3+3*ii+1] );
391  (*rmax_p).push_back( solParam[3+3*ii+2] );
392  }
393  G4double phiTotal = solParam[1];
394  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
395  solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi
396  nplanes, // sections
397  &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0]));
398  }
399  else
400  {
401  std::vector<G4double>* R_c = new std::vector<G4double>;
402  std::vector<G4double>* Z_c = new std::vector<G4double>;
403  for( size_t ii = 0; ii < nplanes; ii++ )
404  {
405  (*R_c).push_back( solParam[3+2*ii] );
406  (*Z_c).push_back( solParam[3+2*ii+1] );
407  }
408  G4double phiTotal = solParam[1];
409  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
410  solid = new G4Polycone( sname, solParam[0], phiTotal, // start,delta-phi
411  nplanes, // sections
412  &((*R_c)[0]), &((*Z_c)[0]));
413  }
414 
415  }
416  else if( stype == "POLYHEDRA" )
417  {
418  size_t nplanes = size_t(solParam[3]);
419  G4bool genericPoly = false;
420  if( solParam.size() == 4+nplanes*3 )
421  {
422  genericPoly = true;
423  }
424  else if( solParam.size() == 4+nplanes*2 )
425  {
426  genericPoly = false;
427  }
428  else
429  {
430  G4String Err1 = "Solid type " + stype + " should have ";
431  G4String Err2 = G4UIcommand::ConvertToString(G4int(4+nplanes*3))
432  + " (Z,Rmin,Rmax)\n";
433  G4String Err3 = "or " + G4UIcommand::ConvertToString(G4int(4+nplanes*2));
434  G4String Err4 = " (RZ corners) parameters,\n";
435  G4String Err5 = "and it has "
436  + G4UIcommand::ConvertToString(G4int(solParam.size()));
437  G4String ErrMessage = Err1 + Err2 + Err3 + Err4 + Err5 + " !";
438  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
439  "InvalidSetup", FatalException, ErrMessage);
440  return 0;
441  }
442 
443  if( genericPoly )
444  {
445  std::vector<G4double>* z_p = new std::vector<G4double>;
446  std::vector<G4double>* rmin_p = new std::vector<G4double>;
447  std::vector<G4double>* rmax_p = new std::vector<G4double>;
448  for( size_t ii = 0; ii < nplanes; ii++ )
449  {
450  (*z_p).push_back( solParam[4+3*ii] );
451  (*rmin_p).push_back( solParam[4+3*ii+1] );
452  (*rmax_p).push_back( solParam[4+3*ii+2] );
453  }
454  G4double phiTotal = solParam[1];
455  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
456  solid = new G4Polyhedra( sname, solParam[0], phiTotal,
457  G4int(solParam[2]), nplanes,
458  &((*z_p)[0]), &((*rmin_p)[0]), &((*rmax_p)[0]));
459  }
460  else
461  {
462  std::vector<G4double>* R_c = new std::vector<G4double>;
463  std::vector<G4double>* Z_c = new std::vector<G4double>;
464  for( size_t ii = 0; ii < nplanes; ii++ )
465  {
466  (*R_c).push_back( solParam[4+2*ii] );
467  (*Z_c).push_back( solParam[4+2*ii+1] );
468  }
469  G4double phiTotal = solParam[1];
470  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
471  solid = new G4Polyhedra( sname, solParam[0], phiTotal,
472  G4int(solParam[2]), nplanes,
473  &((*R_c)[0]), &((*Z_c)[0]));
474  }
475  }
476  else if( stype == "ELLIPTICALTUBE" )
477  {
478  CheckNoSolidParams( stype, 3, solParam.size() );
479  solid = new G4EllipticalTube( sname, solParam[0], solParam[1], solParam[2]);
480  }
481  else if( stype == "ELLIPSOID" )
482  {
483  CheckNoSolidParams( stype, 5, solParam.size() );
484  solid = new G4Ellipsoid( sname, solParam[0], solParam[1], solParam[2],
485  solParam[3], solParam[4] );
486  }
487  else if( stype == "ELLIPTICALCONE" )
488  {
489  CheckNoSolidParams( stype, 4, solParam.size() );
490  solid = new G4EllipticalCone( sname, solParam[0], solParam[1],
491  solParam[2], solParam[3] );
492  }
493  else if( stype == "HYPE" )
494  {
495  CheckNoSolidParams( stype, 5, solParam.size() );
496  solid = new G4Hype( sname, solParam[0], solParam[1], solParam[2],
497  solParam[3], solParam[4] );
498  }
499  else if( stype == "TET" )
500  {
501  CheckNoSolidParams( stype, 12, solParam.size() );
502  G4ThreeVector anchor(solParam[0], solParam[1], solParam[2]);
503  G4ThreeVector p2(solParam[3], solParam[4], solParam[5]);
504  G4ThreeVector p3(solParam[6], solParam[7], solParam[8]);
505  G4ThreeVector p4(solParam[9], solParam[10], solParam[11]);
506  solid = new G4Tet( sname, anchor, p2, p3, p4 );
507  }
508  else if( stype == "TWISTEDBOX" )
509  {
510  CheckNoSolidParams( stype, 4, solParam.size() );
511  solid = new G4TwistedBox( sname, solParam[0], solParam[1],
512  solParam[2], solParam[3]);
513  }
514  else if( stype == "TWISTEDTRAP" )
515  {
516  CheckNoSolidParams( stype, 11, solParam.size() );
517  solid = new G4TwistedTrap( sname, solParam[0], solParam[1], solParam[2],
518  solParam[3], solParam[4], solParam[5], solParam[6],
519  solParam[7], solParam[8], solParam[9], solParam[10] );
520  }
521  else if( stype == "TWISTEDTRD" )
522  {
523  CheckNoSolidParams( stype, 6, solParam.size() );
524  solid = new G4TwistedTrd( sname, solParam[0], solParam[1], solParam[2],
525  solParam[3], solParam[4], solParam[5]);
526  }
527  else if( stype == "TWISTEDTUBS" )
528  {
529  CheckNoSolidParams( stype, 5, solParam.size() );
530  G4double phiTotal = solParam[4];
531  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
532  solid = new G4TwistedTubs( sname, solParam[0], solParam[1], solParam[2],
533  solParam[3], phiTotal);
534  }
535  else if( stype == "BREPBOX" ) // EntityType is = "Closed_Shell"
536  {
537  CheckNoSolidParams( stype, 24, solParam.size() );
538  std::vector<G4Point3D> points;
539  for( size_t ii = 0; ii < 8; ii++ )
540  {
541  points.push_back( G4Point3D(solParam[ii*3+0],
542  solParam[ii*3+1],
543  solParam[ii*3+2]) );
544  }
545  solid = new G4BREPSolidBox( sname, points[0], points[1], points[2],
546  points[3], points[4], points[5], points[6],
547  points[7] );
548  }
549  else if( stype == "BREPCYLINDER" ) // EntityType is = "Closed_Shell"
550  {
551  CheckNoSolidParams( stype, 11, solParam.size() );
552  solid = new G4BREPSolidCylinder( sname,
553  G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
554  G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
555  G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
556  solParam[9], solParam[10] );
557  }
558  else if( stype == "BREPCONE" ) // EntityType is = "Closed_Shell"
559  {
560  CheckNoSolidParams( stype, 12, solParam.size() );
561  solid = new G4BREPSolidCone( sname,
562  G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
563  G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
564  G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
565  solParam[9], solParam[10], solParam[11] );
566  }
567  else if( stype == "BREPSPHERE" ) // EntityType is = "Closed_Shell"
568  {
569  CheckNoSolidParams( stype, 10, solParam.size() );
570  solid = new G4BREPSolidSphere( sname,
571  G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
572  G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
573  G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
574  solParam[9] );
575 
576  }
577  else if( stype == "BREPTORUS" ) // EntityType is = "Closed_Shell"
578  {
579  CheckNoSolidParams( stype, 11, solParam.size() );
580  solid = new G4BREPSolidTorus( sname,
581  G4ThreeVector( solParam[0], solParam[1], solParam[2] ),
582  G4ThreeVector( solParam[3], solParam[4], solParam[5] ),
583  G4ThreeVector( solParam[6], solParam[7], solParam[8] ),
584  solParam[9], solParam[10] );
585  }
586  else if( stype == "BREPPCONE" ) // EntityType is = "Closed_Shell"
587  {
588  size_t nplanes = size_t(solParam[2]);
589  CheckNoSolidParams( stype, 4+3*nplanes, solParam.size() );
590  std::vector<G4double>* z_p = new std::vector<G4double>;
591  std::vector<G4double>* rmin_p = new std::vector<G4double>;
592  std::vector<G4double>* rmax_p = new std::vector<G4double>;
593  for( size_t ii = 0; ii < nplanes; ii++ )
594  {
595  (*z_p).push_back( solParam[4+3*ii] );
596  (*rmin_p).push_back( solParam[4+3*ii+1] );
597  (*rmax_p).push_back( solParam[4+3*ii+2] );
598  }
599  G4double phiTotal = solParam[1];
600  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
601  CheckNoSolidParams( stype, 12, solParam.size() );
602  solid = new G4BREPSolidPCone( sname, solParam[0], phiTotal, // start,dph
603  nplanes, // sections
604  solParam[3], // z_start
605  &((*z_p)[0]), &((*rmin_p)[0]),
606  &((*rmax_p)[0]));
607  }
608  else if( stype == "BREPPOLYHEDRA" ) // EntityType is = "Closed_Shell"
609  {
610  size_t nplanes = size_t(solParam[3]);
611  CheckNoSolidParams( stype, 5+3*nplanes, solParam.size() );
612  std::vector<G4double>* z_p = new std::vector<G4double>;
613  std::vector<G4double>* rmin_p = new std::vector<G4double>;
614  std::vector<G4double>* rmax_p = new std::vector<G4double>;
615  for( size_t ii = 0; ii < nplanes; ii++ )
616  {
617  (*z_p).push_back( solParam[5+3*ii] );
618  (*rmin_p).push_back( solParam[5+3*ii+1] );
619  (*rmax_p).push_back( solParam[5+3*ii+2] );
620  }
621  G4double phiTotal = solParam[1];
622  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
623  CheckNoSolidParams( stype, 12, solParam.size() );
624  solid = new G4BREPSolidPolyhedra( sname, solParam[0], phiTotal, // start,dph
625  G4int(solParam[2]), // sides
626  nplanes, // sections
627  solParam[4], // z_start
628  &((*z_p)[0]), &((*rmin_p)[0]),
629  &((*rmax_p)[0]));
630  }
631  else if( stype == "BREPOPENPCONE" ) // EntityType is = "Closed_Shell"
632  {
633  size_t nplanes = size_t(solParam[2]);
634  std::vector<G4double>* z_p = new std::vector<G4double>;
635  std::vector<G4double>* rmin_p = new std::vector<G4double>;
636  std::vector<G4double>* rmax_p = new std::vector<G4double>;
637  for( size_t ii = 0; ii < nplanes; ii++ )
638  {
639  (*z_p).push_back( solParam[4+3*ii] );
640  (*rmin_p).push_back( solParam[4+3*ii+1] );
641  (*rmax_p).push_back( solParam[4+3*ii+2] );
642  }
643  G4double phiTotal = solParam[1];
644  if( std::fabs(phiTotal - twopi) < angularTolerance ) { phiTotal = twopi; }
645  CheckNoSolidParams( stype, 12, solParam.size() );
646  solid = new G4BREPSolidOpenPCone( sname, solParam[0], phiTotal, // start,dph
647  nplanes, // sections
648  solParam[3], // z_start
649  &((*z_p)[0]), &((*rmin_p)[0]),
650  &((*rmax_p)[0]));
651  }
652  else if( stype == "TESSELLATED" )
653  {
654  G4int nFacets = G4int(solParam[0]);
655  G4int jj = 0;
656  solid = new G4TessellatedSolid(sname);
657  G4TessellatedSolid* solidTS = (G4TessellatedSolid*)(solid);
658  G4VFacet* facet=0;
659 
660  for( G4int ii = 0; ii < nFacets; ii++){
661  G4int nPoints = G4int(solParam[jj+1]);
662  if( G4int(solParam.size()) < jj + nPoints*3 + 2 ) {
663  G4String Err1 = "Too small number of parameters in tesselated solid, it should be at least " + G4UIcommand::ConvertToString(jj + nPoints*3 + 2 );
664  G4String Err2 = " facet number " + G4UIcommand::ConvertToString(ii);
665  G4String Err3 = " number of parameters is " + G4UIcommand::ConvertToString(G4int(solParam.size()));
666  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
667  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
668  "InvalidSetup", FatalException, ErrMessage);
669  return 0;
670  }
671 
672  if( nPoints == 3 )
673  {
674  G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]);
675  G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]);
676  G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]);
677  G4FacetVertexType vertexType = ABSOLUTE;
678  if( solParam[jj+11] == 0 )
679  {
680  vertexType = ABSOLUTE;
681  }
682  else if( solParam[jj+11] == 1 )
683  {
684  vertexType = RELATIVE;
685  }
686  else
687  {
688  G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)";
689  G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
690  G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+11]);
691  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
692  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
693  "InvalidSetup", FatalException, ErrMessage);
694  return 0;
695  }
696  facet = new G4TriangularFacet( pt0, vt1, vt2, vertexType );
697  }
698  else if( nPoints == 4 )
699  {
700  G4ThreeVector pt0(solParam[jj+2],solParam[jj+3],solParam[jj+4]);
701  G4ThreeVector vt1(solParam[jj+5],solParam[jj+6],solParam[jj+7]);
702  G4ThreeVector vt2(solParam[jj+8],solParam[jj+9],solParam[jj+10]);
703  G4ThreeVector vt3(solParam[jj+11],solParam[jj+12],solParam[jj+13]);
704  G4FacetVertexType vertexType = ABSOLUTE;
705  if( solParam[jj+14] == 0 )
706  {
707  vertexType = ABSOLUTE;
708  }
709  else if( solParam[jj+14] == 1 )
710  {
711  vertexType = RELATIVE;
712  }
713  else
714  {
715  G4String Err1 = "Wrong number of vertex type in tesselated solid, it should be 0 =ABSOLUTE) or 1 (=RELATIVE)";
716  G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
717  G4String Err3 = " vertex type is " + G4UIcommand::ConvertToString(solParam[jj+14]);
718  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
719  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
720  "InvalidSetup", FatalException, ErrMessage);
721  return 0;
722  }
723  facet = new G4QuadrangularFacet( pt0, vt1, vt2, vt3, vertexType );
724  }
725  else
726  {
727  G4String Err1 = "Wrong number of points in tesselated solid, it should be 3 or 4";
728  G4String Err2 = " facet number " + G4UIcommand::ConvertToString(G4int(ii));
729  G4String Err3 = " number of points is " + G4UIcommand::ConvertToString(G4int(nPoints));
730  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
731  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
732  "InvalidSetup", FatalException, ErrMessage);
733  return 0;
734  }
735 
736  solidTS->AddFacet( facet );
737  jj += nPoints*3 + 2;
738  }
739 
740  }
741  else if( stype == "EXTRUDED" )
742  {
743  std::vector<G4TwoVector> polygonList;
744  std::vector<G4ExtrudedSolid::ZSection> zsectionList;
745  G4int nPolygons = G4int(solParam[0]);
746  G4int ii = 1;
747  G4int nMax = nPolygons*2+1;
748  for( ;ii < nMax; ii+=2 )
749  {
750  polygonList.push_back( G4TwoVector(solParam[ii],solParam[ii+1]) );
751  }
752  G4int nZSections = G4int(solParam[ii]);
753  nMax = nPolygons*2 + nZSections*4 + 2;
754  ii++;
755  for( ; ii < nMax; ii+=4 )
756  {
757  G4TwoVector offset(solParam[ii+1],solParam[ii+2]);
758  zsectionList.push_back( G4ExtrudedSolid::ZSection(solParam[ii],offset,solParam[ii+3]) );
759  }
760  solid = new G4ExtrudedSolid( sname, polygonList, zsectionList );
761 
762  }
763  else if( stype.substr(0,7) == "Boolean" )
764  {
765  const G4tgrSolidBoolean* solb = dynamic_cast<const G4tgrSolidBoolean*>(sol);
766  if (!solb)
767  {
768  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
769  "InvalidSetup", FatalException, "Invalid Solid pointer");
770  return 0;
771  }
772  G4VSolid* sol1 = FindOrConstructG4Solid( solb->GetSolid(0));
773  G4VSolid* sol2 = FindOrConstructG4Solid( solb->GetSolid(1));
776  G4ThreeVector relPlace = solb->GetRelativePlace();
777 
778  if( stype == "Boolean_UNION" )
779  {
780  solid = new G4UnionSolid( sname, sol1, sol2, relRotMat, relPlace );
781  }
782  else if( stype == "Boolean_SUBTRACTION" )
783  {
784  solid = new G4SubtractionSolid( sname, sol1, sol2, relRotMat, relPlace );
785  }
786  else if( stype == "Boolean_INTERSECTION" )
787  {
788  solid = new G4IntersectionSolid( sname, sol1, sol2, relRotMat, relPlace );
789  }
790  else
791  {
792  G4String ErrMessage = "Unknown Boolean type " + stype;
793  G4Exception("G4tgbVolume::FindOrConstructG4Solid()",
794  "InvalidSetup", FatalException, ErrMessage);
795  return 0;
796  }
797  }
798  else
799  {
800  G4String ErrMessage = "Solids of type " + stype
801  + " not implemented yet, sorry...";
802  G4Exception("G4tgbVolume::FindOrConstructG4Solid()", "NotImplemented",
803  FatalException, ErrMessage);
804  return 0;
805  }
806 
807 #ifdef G4VERBOSE
809  {
810  G4cout << " G4tgbVolume::FindOrConstructG4Solid()" << G4endl
811  << " Created solid " << sname
812  << " of type " << solid->GetEntityType() << G4endl;
813  }
814 #endif
815 
816 #ifdef G4VERBOSE
818  {
819  G4cout << " Constructing new G4Solid: "
820  << *solid << G4endl;
821  }
822 #endif
823 
824  return solid;
825 }
826 
827 //-------------------------------------------------------------------
829  const unsigned int NoParamExpected,
830  const unsigned int NoParam )
831 {
832  if( NoParamExpected != NoParam )
833  {
834  G4String Err1 = "Solid type " + solidType + " should have ";
835  G4String Err2 = G4UIcommand::ConvertToString(G4int(NoParamExpected))
836  + " parameters,\n";
837  G4String Err3 = "and it has "
839  G4String ErrMessage = Err1 + Err2 + Err3 + " !";
840  G4Exception("G4tgbVolume::CheckNoSolidParams()", "InvalidSetup",
841  FatalException, ErrMessage);
842  }
843 }
844 
845 
846 //-------------------------------------------------------------------
848 {
849  G4LogicalVolume* logvol;
850 
851 #ifdef G4VERBOSE
853  {
854  G4cout << " G4tgbVolume::ConstructG4LogVol() - " << GetName() << G4endl;
855  }
856 #endif
857 
858  //----------- Get the material first
860  ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
861  if( mate == 0 )
862  {
863  G4String ErrMessage = "Material not found "
864  + theTgrVolume->GetMaterialName()
865  + " for volume " + GetName() + ".";
866  G4Exception("G4tgbVolume::ConstructG4LogVol()", "InvalidSetup",
867  FatalException, ErrMessage);
868  }
869 #ifdef G4VERBOSE
871  {
872  G4cout << " G4tgbVolume::ConstructG4LogVol() -"
873  << " Material constructed: " << mate->GetName() << G4endl;
874  }
875 #endif
876 
877  //---------- Construct the LV
878  logvol = new G4LogicalVolume( const_cast<G4VSolid*>(solid),
879  const_cast<G4Material*>(mate), GetName() );
880 
881 #ifdef G4VERBOSE
883  {
884  G4cout << " Constructing new G4LogicalVolume: "
885  << logvol->GetName() << " mate " << mate->GetName() << G4endl;
886  }
887 #endif
888 
889  //---------- Set visibility and colour
890  if( !GetVisibility() || GetColour()[0] != -1 )
891  {
892  G4VisAttributes* visAtt = new G4VisAttributes();
893 #ifdef G4VERBOSE
895  {
896  G4cout << " Constructing new G4VisAttributes: "
897  << *visAtt << G4endl;
898  }
899 #endif
900 
901  if( !GetVisibility() )
902  {
903  visAtt->SetVisibility( GetVisibility() );
904  }
905  else if( GetColour()[0] != -1 )
906  {
907  // this else should not be necessary, because if the visibility
908  // is set to off, colour should have no effect. But it does not
909  // work: if you set colour and vis off, it is visualized!?!?!?
910 
911  const G4double* col = GetColour();
912  if( col[3] == -1. )
913  {
914  visAtt->SetColour( G4Colour(col[0],col[1],col[2]));
915  }
916  else
917  {
918  visAtt->SetColour( G4Colour(col[0],col[1],col[2],col[3]));
919  }
920  }
921  logvol->SetVisAttributes(visAtt);
922  }
923 
924 #ifdef G4VERBOSE
926  {
927  G4cout << " G4tgbVolume::ConstructG4LogVol() -"
928  << " Created logical volume: " << GetName() << G4endl;
929  }
930 #endif
931 
932  return logvol;
933 }
934 
935 
936 //-------------------------------------------------------------------
939  const G4LogicalVolume* currentLV,
940  const G4LogicalVolume* parentLV )
941 {
942  G4VPhysicalVolume* physvol = 0;
943  G4int copyNo;
944 
945  //----- Case of placement of top volume
946  if( place == 0 )
947  {
948 #ifdef G4VERBOSE
950  {
951  G4cout << " G4tgbVolume::ConstructG4PhysVol() - World: "
952  << GetName() << G4endl;
953  }
954 #endif
955  physvol = new G4PVPlacement(0, G4ThreeVector(),
956  const_cast<G4LogicalVolume*>(currentLV),
957  GetName(), 0, false, 0,
958  theTgrVolume->GetCheckOverlaps());
959 #ifdef G4VERBOSE
961  {
962  G4cout << " Constructing new : G4PVPlacement "
963  << physvol->GetName() << G4endl;
964  }
965 #endif
966  }
967  else
968  {
969  copyNo = place->GetCopyNo();
970 
971 #ifdef G4VERBOSE
973  {
974  G4cout << " G4tgbVolume::ConstructG4PhysVol() - " << GetName() << G4endl
975  << " inside " << parentLV->GetName() << " copy No: " << copyNo
976  << " of type= " << theTgrVolume->GetType() << G4endl
977  << " placement type= " << place->GetType() << G4endl;
978  }
979 #endif
980 
981  if( theTgrVolume->GetType() == "VOLSimple" )
982  {
983  //----- Get placement
984 #ifdef G4VERBOSE
986  {
987  G4cout << " G4tgbVolume::ConstructG4PhysVol() - Placement type = "
988  << place->GetType() << G4endl;
989  }
990 #endif
991 
992  //--------------- If it is G4tgrPlaceSimple
993  if( place->GetType() == "PlaceSimple" )
994  {
995  //----- Get rotation matrix
996  G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
997  G4String rmName = placeSimple->GetRotMatName();
998 
1000  ->FindOrBuildG4RotMatrix( rmName );
1001  //----- Place volume in mother
1002  G4double check = (rotmat->colX().cross(rotmat->colY()))*rotmat->colZ();
1003  G4double tol = 1.0e-3;
1004  //---- Check that matrix is ortogonal
1005  if (1-std::abs(check)>tol)
1006  {
1007  G4cerr << " Matrix : " << rmName << " " << rotmat->colX()
1008  << " " << rotmat->colY() << " " << rotmat->colZ() << G4endl
1009  << " product x X y * z = " << check << " x X y "
1010  << rotmat->colX().cross(rotmat->colY()) << G4endl;
1011  G4String ErrMessage = "Rotation is not ortogonal " + rmName + " !";
1012  G4Exception("G4tgbVolume::ConstructG4PhysVol()",
1013  "InvalidSetup", FatalException, ErrMessage);
1014  //---- Check if it is reflection
1015  }
1016  else if (1+check<=tol)
1017  {
1018  G4Translate3D transl = place->GetPlacement();
1019  G4Transform3D trfrm = transl * G4Rotate3D(*rotmat);
1020  physvol = (G4ReflectionFactory::Instance()->Place(trfrm, GetName(),
1021  const_cast<G4LogicalVolume*>(currentLV),
1022  const_cast<G4LogicalVolume*>(parentLV),
1023  false, copyNo, false )).first;
1024  }
1025  else
1026  {
1027 #ifdef G4VERBOSE
1028  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1029  {
1030  G4cout << "Construction new G4VPhysicalVolume"
1031  << " through G4ReflectionFactory " << GetName()
1032  << " in volume " << parentLV->GetName()
1033  << " copyNo " << copyNo
1034  << " position " << place->GetPlacement()
1035  << " ROT " << rotmat->colX()
1036  << " " << rotmat->colY()
1037  << " " << rotmat->colZ() << G4endl;
1038  }
1039 #endif
1040  physvol = new G4PVPlacement( rotmat, place->GetPlacement(),
1041  const_cast<G4LogicalVolume*>(currentLV),
1042  GetName(),
1043  const_cast<G4LogicalVolume*>(parentLV),
1044  false, copyNo,
1045  theTgrVolume->GetCheckOverlaps());
1046  }
1047 
1048  //--------------- If it is G4tgrPlaceParam
1049  }
1050  else if( place->GetType() == "PlaceParam" )
1051  {
1053 
1054  //----- See what parameterisation type
1055 #ifdef G4VERBOSE
1056  if( G4tgrMessenger::GetVerboseLevel() >= 2 )
1057  {
1058  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1059  << " param: " << GetName() << " in " << parentLV->GetName()
1060  << " param type= " << dp->GetParamType() << G4endl;
1061  }
1062 #endif
1063 
1064  G4tgbPlaceParameterisation * param=0;
1065 
1066  if( (dp->GetParamType() == "CIRCLE")
1067  || (dp->GetParamType() == "CIRCLE_XY")
1068  || (dp->GetParamType() == "CIRCLE_XZ")
1069  || (dp->GetParamType() == "CIRCLE_YZ") )
1070  {
1071  param = new G4tgbPlaceParamCircle(dp);
1072 
1073  }
1074  else if( (dp->GetParamType() == "LINEAR")
1075  || (dp->GetParamType() == "LINEAR_X")
1076  || (dp->GetParamType() == "LINEAR_Y")
1077  || (dp->GetParamType() == "LINEAR_Z") )
1078  {
1079  param = new G4tgbPlaceParamLinear(dp);
1080 
1081  }
1082  else if( (dp->GetParamType() == "SQUARE")
1083  || (dp->GetParamType() == "SQUARE_XY")
1084  || (dp->GetParamType() == "SQUARE_XZ")
1085  || (dp->GetParamType() == "SQUARE_YZ") )
1086  {
1087  param = new G4tgbPlaceParamSquare(dp);
1088  }
1089  else
1090  {
1091  G4String ErrMessage = "Parameterisation has wrong type, TYPE: "
1092  + G4String(dp->GetParamType()) + " !";
1093  G4Exception("G4tgbVolume::ConstructG4PhysVol", "WrongArgument",
1094  FatalException, ErrMessage);
1095  return 0;
1096  }
1097 #ifdef G4VERBOSE
1098  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1099  {
1100  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1101  << " New G4PVParameterised: " << GetName() << " vol "
1102  << currentLV->GetName() << " in vol " << parentLV->GetName()
1103  << " axis " << param->GetAxis() << " nCopies "
1104  << param->GetNCopies() << G4endl;
1105  }
1106 #endif
1107  physvol = new G4PVParameterised(GetName(),
1108  const_cast<G4LogicalVolume*>(currentLV),
1109  const_cast<G4LogicalVolume*>(parentLV),
1110  EAxis(param->GetAxis()),
1111  param->GetNCopies(), param);
1112 #ifdef G4VERBOSE
1113  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1114  {
1115  G4cout << " Constructing new G4PVParameterised: "
1116  << physvol->GetName() << " in volume " << parentLV->GetName()
1117  << " N copies " << param->GetNCopies()
1118  << " axis " << param->GetAxis() << G4endl;
1119  }
1120 #endif
1121 
1122  }
1123  else if( place->GetType() == "PlaceReplica" )
1124  {
1125  //--------------- If it is PlaceReplica
1126  G4tgrPlaceDivRep* dpr = (G4tgrPlaceDivRep*)place;
1127 
1128 #ifdef G4VERBOSE
1129  if( G4tgrMessenger::GetVerboseLevel() >= 2 )
1130  {
1131  G4cout << " G4tgbVolume::ConstructG4PhysVol() -" << G4endl
1132  << " replica" << " " << currentLV->GetName()
1133  << " in " << parentLV->GetName()
1134  << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1135  << " offset " << dpr->GetOffset() << G4endl;
1136  }
1137 #endif
1138  physvol = new G4PVReplica(GetName(),
1139  const_cast<G4LogicalVolume*>(currentLV),
1140  const_cast<G4LogicalVolume*>(parentLV),
1141  EAxis(dpr->GetAxis()), dpr->GetNDiv(),
1142  dpr->GetWidth(), dpr->GetOffset());
1143 #ifdef G4VERBOSE
1144  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1145  {
1146  G4cout << " Constructing new G4PVReplica: "
1147  << currentLV->GetName()
1148  << " in " << parentLV->GetName()
1149  << " NDiv " << dpr->GetNDiv() << " Width " << dpr->GetWidth()
1150  << " offset " << dpr->GetOffset() << G4endl;
1151  }
1152 #endif
1153  }
1154  }
1155  else if( theTgrVolume->GetType() == "VOLDivision" )
1156  {
1157  G4tgrVolumeDivision* volr = (G4tgrVolumeDivision*)theTgrVolume;
1158  G4tgrPlaceDivRep* placeDiv = volr->GetPlaceDivision() ;
1159  G4VSolid* solid = BuildSolidForDivision( parentLV->GetSolid(), placeDiv->GetAxis() );
1161  ->FindOrBuildG4Material( theTgrVolume->GetMaterialName() );
1162  G4LogicalVolume* divLV = new G4LogicalVolume(solid,
1163  const_cast<G4Material*>(mate),
1164  GetName() );
1165 #ifdef G4VERBOSE
1166  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1167  {
1168  G4cout << " Constructed new G4LogicalVolume for division: "
1169  << divLV->GetName() << " mate " << mate->GetName() << G4endl;
1170  }
1171 #endif
1172 
1173  G4DivType divType = placeDiv->GetDivType();
1174  switch (divType)
1175  {
1176  case DivByNdiv:
1177  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1178  const_cast<G4LogicalVolume*>(parentLV),
1179  placeDiv->GetAxis(), placeDiv->GetNDiv(),
1180  placeDiv->GetOffset());
1181 #ifdef G4VERBOSE
1182  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1183  {
1184  G4cout << " Constructing new G4PVDivision by number of divisions: "
1185  << GetName() << " in " << parentLV->GetName()
1186  << " axis " << placeDiv->GetAxis()
1187  << " Ndiv " << placeDiv->GetNDiv()
1188  << " offset " << placeDiv->GetOffset() << G4endl;
1189  }
1190 #endif
1191  break;
1192  case DivByWidth:
1193  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1194  const_cast<G4LogicalVolume*>(parentLV),
1195  placeDiv->GetAxis(), placeDiv->GetWidth(),
1196  placeDiv->GetOffset());
1197 #ifdef G4VERBOSE
1198  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1199  {
1200  G4cout << " Constructing new G4PVDivision by width: "
1201  << GetName() << " in " << parentLV->GetName()
1202  << " axis " << placeDiv->GetAxis()
1203  << " width " << placeDiv->GetWidth()
1204  << " offset " << placeDiv->GetOffset() << G4endl;
1205  }
1206 #endif
1207  break;
1208  case DivByNdivAndWidth:
1209  physvol = new G4PVDivision(GetName(), (G4LogicalVolume*)divLV,
1210  const_cast<G4LogicalVolume*>(parentLV),
1211  placeDiv->GetAxis(), placeDiv->GetNDiv(),
1212  placeDiv->GetWidth(),
1213  placeDiv->GetOffset());
1214 #ifdef G4VERBOSE
1215  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1216  {
1217  G4cout << " Constructing new G4PVDivision by width"
1218  << " and number of divisions: "
1219  << GetName() << " in " << parentLV->GetName()
1220  << " axis " << placeDiv->GetAxis()
1221  << " Ndiv " << placeDiv->GetNDiv()
1222  << " width " << placeDiv->GetWidth()
1223  << " offset " << placeDiv->GetOffset() << G4endl;
1224  }
1225 #endif
1226  break;
1227  }
1228  }
1229  else if( theTgrVolume->GetType() == "VOLAssembly" )
1230  {
1231  // Define one layer as one assembly volume
1232  G4tgrVolumeAssembly * tgrAssembly = (G4tgrVolumeAssembly *)theTgrVolume;
1233 
1234  if( !theG4AssemblyVolume )
1235  {
1236  theG4AssemblyVolume = new G4AssemblyVolume;
1237 #ifdef G4VERBOSE
1238  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1239  {
1240  G4cout << " Constructing new G4AssemblyVolume: "
1241  << " number of assembly components "
1242  << tgrAssembly->GetNoComponents() << G4endl;
1243  }
1244 #endif
1246  for( G4int ii = 0; ii < tgrAssembly->GetNoComponents(); ii++ )
1247  {
1248  // Rotation and translation of a plate inside the assembly
1249 
1250  G4ThreeVector transl = tgrAssembly->GetComponentPos(ii);
1251  G4String rmName = tgrAssembly->GetComponentRM(ii);
1253  ->FindOrBuildG4RotMatrix( rmName );
1254 
1255  //----- Get G4LogicalVolume of component
1256  G4String lvname = tgrAssembly->GetComponentName(ii);
1257  G4LogicalVolume* logvol = g4vmgr->FindG4LogVol( lvname);
1258  if( logvol == 0 )
1259  {
1260  g4vmgr->FindVolume( lvname )->ConstructG4Volumes( 0, 0);
1261  logvol = g4vmgr->FindG4LogVol( lvname, true );
1262  }
1263  // Fill the assembly by the plates
1264  theG4AssemblyVolume->AddPlacedVolume( logvol, transl, rotmat );
1265 #ifdef G4VERBOSE
1266  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1267  {
1268  G4cout << " G4AssemblyVolume->AddPlacedVolume " << ii
1269  << " " << logvol->GetName()
1270  << " translation " << transl
1271  << " rotmat " << rotmat->colX()
1272  << " " << rotmat->colY()
1273  << " " << rotmat->colZ() << G4endl;
1274  }
1275 #endif
1276  }
1277  }
1278 
1279  // Rotation and Translation of the assembly inside the world
1280 
1281  G4tgrPlaceSimple* placeSimple = (G4tgrPlaceSimple*)place;
1282  G4String rmName = placeSimple->GetRotMatName();
1284  ->FindOrBuildG4RotMatrix( rmName );
1285  G4ThreeVector transl = place->GetPlacement();
1286 
1287  G4LogicalVolume* parentLV_nonconst =
1288  const_cast<G4LogicalVolume*>(parentLV);
1289  theG4AssemblyVolume->MakeImprint( parentLV_nonconst, transl, rotmat );
1290 
1291  }
1292  else // If it is G4tgrVolumeAssembly
1293  {
1294  G4String ErrMessage = "Volume type not supported: "
1295  + theTgrVolume->GetType() + ", sorry...";
1296  G4Exception("G4tgbVolume::ConstructG4PhysVol()", "NotImplemented",
1297  FatalException, ErrMessage);
1298  }
1299  }
1300 
1301  return physvol;
1302 }
1303 
1304 
1305 //-------------------------------------------------------------------
1307 {
1308  G4VSolid* solid=0;
1309  G4double redf = (parentSolid->GetExtent().GetXmax()-parentSolid->GetExtent().GetXmin());
1310  redf = std::min(redf,parentSolid->GetExtent().GetYmax()-parentSolid->GetExtent().GetYmin());
1311  redf = std::min(redf,parentSolid->GetExtent().GetZmax()-parentSolid->GetExtent().GetZmin());
1312  redf *= 0.001; //make daugther much smaller, to fit in parent
1313 
1314  if( parentSolid->GetEntityType() == "G4Box" )
1315  {
1316  G4Box* psolid = (G4Box*)(parentSolid);
1317  solid = new G4Box(GetName(), psolid->GetXHalfLength()*redf,
1318  psolid->GetZHalfLength()*redf,
1319  psolid->GetZHalfLength()*redf);
1320  }
1321  else if ( parentSolid->GetEntityType() == "G4Tubs" )
1322  {
1323  G4Tubs* psolid = (G4Tubs*)(parentSolid);
1324  solid = new G4Tubs( GetName(), psolid->GetInnerRadius()*redf,
1325  psolid->GetOuterRadius()*redf,
1326  psolid->GetZHalfLength()*redf,
1327  psolid->GetSPhi(), psolid->GetDPhi());
1328  }
1329  else if ( parentSolid->GetEntityType() == "G4Cons" )
1330  {
1331  G4Cons* psolid = (G4Cons*)(parentSolid);
1332  solid = new G4Cons( GetName(), psolid->GetInnerRadiusMinusZ()*redf,
1333  psolid->GetOuterRadiusMinusZ()*redf,
1334  psolid->GetInnerRadiusPlusZ()*redf,
1335  psolid->GetOuterRadiusPlusZ()*redf,
1336  psolid->GetZHalfLength()*redf,
1337  psolid->GetSPhi(), psolid->GetDPhi());
1338  }
1339  else if ( parentSolid->GetEntityType() == "G4Trd" )
1340  {
1341  G4Trd* psolid = (G4Trd*)(parentSolid);
1342  G4double mpDx1 = psolid->GetXHalfLength1();
1343  G4double mpDx2 = psolid->GetXHalfLength2();
1344 
1345  if( axis == kXAxis && std::fabs(mpDx1 - mpDx2) > G4GeometryTolerance::GetInstance()->GetSurfaceTolerance() )
1346  {
1347  solid = new G4Trap( GetName(), psolid->GetZHalfLength()*redf,
1348  psolid->GetYHalfLength1()*redf,
1349  psolid->GetXHalfLength2()*redf,
1350  psolid->GetXHalfLength1()*redf );
1351  }
1352  else
1353  {
1354  solid = new G4Trd( GetName(), psolid->GetXHalfLength1()*redf,
1355  psolid->GetXHalfLength2()*redf,
1356  psolid->GetYHalfLength1()*redf,
1357  psolid->GetYHalfLength2()*redf,
1358  psolid->GetZHalfLength()*redf);
1359  }
1360 
1361  }
1362  else if ( parentSolid->GetEntityType() == "G4Para" )
1363  {
1364  G4Para* psolid = (G4Para*)(parentSolid);
1365  solid = new G4Para( GetName(), psolid->GetXHalfLength()*redf,
1366  psolid->GetYHalfLength()*redf,
1367  psolid->GetZHalfLength()*redf,
1368  std::atan(psolid->GetTanAlpha()),
1369  psolid->GetSymAxis().theta(),
1370  psolid->GetSymAxis().phi() );
1371  }
1372  else if ( parentSolid->GetEntityType() == "G4Polycone" )
1373  {
1374  G4Polycone* psolid = (G4Polycone*)(parentSolid);
1375  G4PolyconeHistorical origParam = *(psolid->GetOriginalParameters());
1376  for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1377  {
1378  origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1379  origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1380  }
1381  solid = new G4Polycone( GetName(), psolid->GetStartPhi(),
1382  psolid->GetEndPhi(),
1383  origParam.Num_z_planes, origParam.Z_values,
1384  origParam.Rmin, origParam.Rmax);
1385 
1386  }
1387  else if ( parentSolid->GetEntityType() == "G4Polyhedra" )
1388  {
1389  G4Polyhedra* psolid = (G4Polyhedra*)(parentSolid);
1390  G4PolyhedraHistorical origParam = *(psolid->GetOriginalParameters());
1391  for( G4int ii = 0; ii < origParam.Num_z_planes; ii++ )
1392  {
1393  origParam.Rmin[ii] = origParam.Rmin[ii]*redf;
1394  origParam.Rmax[ii] = origParam.Rmax[ii]*redf;
1395  }
1396  solid = new G4Polyhedra( GetName(), psolid->GetStartPhi(),
1397  psolid->GetEndPhi(),
1398  psolid->GetNumSide(),
1399  origParam.Num_z_planes, origParam.Z_values,
1400  origParam.Rmin, origParam.Rmax);
1401  }
1402  else
1403  {
1404  G4String ErrMessage = "Solid type not supported. VOLUME= " + GetName()
1405  + " Solid type= " + parentSolid->GetEntityType() + "\n"
1406  + "Only supported types are: G4Box, G4Tubs, G4Cons,"
1407  + " G4Trd, G4Para, G4Polycone, G4Polyhedra.";
1408  G4Exception("G4tgbVolume::BuildSolidForDivision()", "NotImplemented",
1409  FatalException, ErrMessage);
1410  return 0;
1411  }
1412 
1413 #ifdef G4VERBOSE
1414  if( G4tgrMessenger::GetVerboseLevel() >= 1 )
1415  {
1416  G4cout << " Constructing new G4Solid for division: "
1417  << *solid << G4endl;
1418  }
1419 #endif
1420  return solid;
1421 }
1422