source: trunk/source/persistency/ascii/src/G4tgbVolume.cc@ 1127

Last change on this file since 1127 was 1035, checked in by garnier, 17 years ago

dossiers oublies

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