source: trunk/source/geometry/divisions/src/G4ParameterisationPolycone.cc@ 1044

Last change on this file since 1044 was 1010, checked in by garnier, 17 years ago

update

File size: 17.9 KB
Line 
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: G4ParameterisationPolycone.cc,v 1.15 2006/06/29 18:18:44 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// class G4ParameterisationPolycone Implementation file
31//
32// 26.05.03 - P.Arce, Initial version
33// 08.04.04 - I.Hrivnacova, Implemented reflection
34//---------------------------------------------------------------------
35
36#include "G4ParameterisationPolycone.hh"
37
38#include <iomanip>
39#include "G4ThreeVector.hh"
40#include "G4RotationMatrix.hh"
41#include "G4VPhysicalVolume.hh"
42#include "G4LogicalVolume.hh"
43#include "G4ReflectedSolid.hh"
44
45//-----------------------------------------------------------------------
46G4VParameterisationPolycone::
47G4VParameterisationPolycone( EAxis axis, G4int nDiv, G4double width,
48 G4double offset, G4VSolid* msolid,
49 DivisionType divType )
50 : G4VDivisionParameterisation( axis, nDiv, width, offset, divType, msolid )
51{
52 G4Polycone* msol = (G4Polycone*)(msolid);
53 if ((msolid->GetEntityType() != "G4ReflectedSolid") && (msol->IsGeneric()))
54 {
55 G4String message =
56 "Sorry, generic construct for G4Polycone NOT supported.\n Solid: "
57 + msol->GetName();
58 G4Exception("G4VParameterisationPolycone::G4VParameterisationPolycone()",
59 "NotSupported", FatalException, message);
60 }
61 if (msolid->GetEntityType() == "G4ReflectedSolid")
62 {
63 // Get constituent solid
64 G4VSolid* mConstituentSolid
65 = ((G4ReflectedSolid*)msolid)->GetConstituentMovedSolid();
66 msol = (G4Polycone*)(mConstituentSolid);
67
68 // Get parameters
69 G4int nofZplanes = msol->GetOriginalParameters()->Num_z_planes;
70 G4double* zValues = msol->GetOriginalParameters()->Z_values;
71 G4double* rminValues = msol->GetOriginalParameters()->Rmin;
72 G4double* rmaxValues = msol->GetOriginalParameters()->Rmax;
73
74 // Invert z values
75 G4double* zValuesRefl = new double[nofZplanes];
76 for (G4int i=0; i<nofZplanes; i++) zValuesRefl[i] = - zValues[i];
77
78 G4Polycone* newSolid
79 = new G4Polycone(msol->GetName(),
80 msol->GetStartPhi(),
81 msol->GetEndPhi() - msol->GetStartPhi(),
82 nofZplanes, zValuesRefl, rminValues, rmaxValues);
83
84 delete [] zValuesRefl;
85
86 msol = newSolid;
87 fmotherSolid = newSolid;
88 fReflectedSolid = true;
89 fDeleteSolid = true;
90 }
91}
92
93//---------------------------------------------------------------------
94G4VParameterisationPolycone::~G4VParameterisationPolycone()
95{
96}
97
98//---------------------------------------------------------------------
99G4ParameterisationPolyconeRho::
100G4ParameterisationPolyconeRho( EAxis axis, G4int nDiv,
101 G4double width, G4double offset,
102 G4VSolid* msolid, DivisionType divType )
103 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
104{
105 CheckParametersValidity();
106 SetType( "DivisionPolyconeRho" );
107
108 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
109 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
110
111 if( divType == DivWIDTH )
112 {
113 fnDiv = CalculateNDiv( origparamMother->Rmax[0]
114 - origparamMother->Rmin[0], width, offset );
115 }
116 else if( divType == DivNDIV )
117 {
118 fwidth = CalculateWidth( origparamMother->Rmax[0]
119 - origparamMother->Rmin[0], nDiv, offset );
120 }
121
122#ifdef G4DIVDEBUG
123 if( verbose >= -1 )
124 {
125 G4cout << " G4ParameterisationPolyconeRho - # divisions " << fnDiv
126 << " = " << nDiv << G4endl
127 << " Offset " << foffset << " = " << offset << G4endl
128 << " Width " << fwidth << " = " << width << G4endl;
129 }
130#endif
131}
132
133//---------------------------------------------------------------------
134G4ParameterisationPolyconeRho::~G4ParameterisationPolyconeRho()
135{
136}
137
138//---------------------------------------------------------------------
139void G4ParameterisationPolyconeRho::CheckParametersValidity()
140{
141 G4VDivisionParameterisation::CheckParametersValidity();
142
143 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
144
145 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
146 {
147 G4cerr << "WARNING - "
148 << "G4ParameterisationPolyconeRho::CheckParametersValidity()"
149 << G4endl
150 << " Solid " << msol->GetName() << G4endl
151 << " Division along R will be done with a width "
152 << "different for each solid section." << G4endl
153 << " WIDTH will not be used !" << G4endl;
154 }
155 if( foffset != 0. )
156 {
157 G4cerr << "WARNING - "
158 << "G4ParameterisationPolyconeRho::CheckParametersValidity()"
159 << G4endl
160 << " Solid " << msol->GetName() << G4endl
161 << " Division along R will be done with a width "
162 << "different for each solid section." << G4endl
163 << " OFFSET will not be used !" << G4endl;
164 }
165
166}
167
168//------------------------------------------------------------------------
169G4double G4ParameterisationPolyconeRho::GetMaxParameter() const
170{
171 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
172 G4PolyconeHistorical* original_pars = msol->GetOriginalParameters();
173 return original_pars->Rmax[0] - original_pars->Rmin[0];
174}
175
176
177//---------------------------------------------------------------------
178void
179G4ParameterisationPolyconeRho::
180ComputeTransformation( const G4int, G4VPhysicalVolume* physVol ) const
181{
182 //----- translation
183 G4ThreeVector origin(0.,0.,0.);
184 //----- set translation
185 physVol->SetTranslation( origin );
186
187 //----- calculate rotation matrix: unit
188
189#ifdef G4DIVDEBUG
190 if( verbose >= 2 )
191 {
192 G4cout << " G4ParameterisationPolyconeRho " << G4endl
193 << " foffset: " << foffset
194 << " - fwidth: " << fwidth << G4endl;
195 }
196#endif
197
198 ChangeRotMatrix( physVol );
199
200#ifdef G4DIVDEBUG
201 if( verbose >= 2 )
202 {
203 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeRho "
204 << G4endl
205 << " Position: " << origin/mm
206 << " - Width: " << fwidth/deg
207 << " - Axis: " << faxis << G4endl;
208 }
209#endif
210}
211
212//---------------------------------------------------------------------
213void
214G4ParameterisationPolyconeRho::
215ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
216 const G4VPhysicalVolume* ) const
217{
218 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
219
220 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
221 G4PolyconeHistorical origparam( *origparamMother );
222 G4int nZplanes = origparamMother->Num_z_planes;
223
224 G4double width = 0.;
225 for( G4int ii = 0; ii < nZplanes; ii++ )
226 {
227 width = CalculateWidth( origparamMother->Rmax[ii]
228 - origparamMother->Rmin[ii], fnDiv, foffset );
229 origparam.Rmin[ii] = origparamMother->Rmin[ii]+foffset+width*copyNo;
230 origparam.Rmax[ii] = origparamMother->Rmin[ii]+foffset+width*(copyNo+1);
231 }
232
233 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
234 pcone.Reset(); // reset to new solid parameters
235
236#ifdef G4DIVDEBUG
237 if( verbose >= -2 )
238 {
239 G4cout << "G4ParameterisationPolyconeRho::ComputeDimensions()" << G4endl
240 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
241 pcone.DumpInfo();
242 }
243#endif
244}
245
246//---------------------------------------------------------------------
247G4ParameterisationPolyconePhi::
248G4ParameterisationPolyconePhi( EAxis axis, G4int nDiv,
249 G4double width, G4double offset,
250 G4VSolid* msolid, DivisionType divType )
251 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
252{
253 CheckParametersValidity();
254 SetType( "DivisionPolyconePhi" );
255
256 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
257 G4double deltaPhi = msol->GetEndPhi() - msol->GetStartPhi();
258
259 if( divType == DivWIDTH )
260 {
261 fnDiv = CalculateNDiv( deltaPhi, width, offset );
262 }
263 else if( divType == DivNDIV )
264 {
265 fwidth = CalculateWidth( deltaPhi, nDiv, offset );
266 }
267
268#ifdef G4DIVDEBUG
269 if( verbose >= 1 )
270 {
271 G4cout << " G4ParameterisationPolyconePhi - # divisions " << fnDiv
272 << " = " << nDiv << G4endl
273 << " Offset " << foffset/deg << " = " << offset/deg << G4endl
274 << " Width " << fwidth/deg << " = " << width/deg << G4endl;
275 }
276#endif
277}
278
279//---------------------------------------------------------------------
280G4ParameterisationPolyconePhi::~G4ParameterisationPolyconePhi()
281{
282}
283
284//------------------------------------------------------------------------
285G4double G4ParameterisationPolyconePhi::GetMaxParameter() const
286{
287 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
288 return msol->GetEndPhi() - msol->GetStartPhi();
289}
290
291//---------------------------------------------------------------------
292void
293G4ParameterisationPolyconePhi::
294ComputeTransformation( const G4int copyNo, G4VPhysicalVolume *physVol ) const
295{
296 //----- translation
297 G4ThreeVector origin(0.,0.,0.);
298 //----- set translation
299 physVol->SetTranslation( origin );
300
301 //----- calculate rotation matrix (so that all volumes point to the centre)
302 G4double posi = foffset + copyNo*fwidth;
303
304#ifdef G4DIVDEBUG
305 if( verbose >= 2 )
306 {
307 G4cout << " G4ParameterisationPolyconePhi - position: " << posi/deg
308 << G4endl
309 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
310 << " - fwidth: " << fwidth/deg << G4endl;
311 }
312#endif
313
314 ChangeRotMatrix( physVol, -posi );
315
316#ifdef G4DIVDEBUG
317 if( verbose >= 2 )
318 {
319 G4cout << std::setprecision(8) << " G4ParameterisationPolyconePhi "
320 << copyNo << G4endl
321 << " Position: " << origin << " - Width: " << fwidth
322 << " - Axis: " << faxis << G4endl;
323 }
324#endif
325}
326
327//---------------------------------------------------------------------
328void
329G4ParameterisationPolyconePhi::
330ComputeDimensions( G4Polycone& pcone, const G4int,
331 const G4VPhysicalVolume* ) const
332{
333 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
334
335 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
336 G4PolyconeHistorical origparam( *origparamMother );
337 origparam.Start_angle = origparamMother->Start_angle;
338 origparam.Opening_angle = fwidth;
339
340 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
341 pcone.Reset(); // reset to new solid parameters
342
343#ifdef G4DIVDEBUG
344 if( verbose >= 2 )
345 {
346 G4cout << "G4ParameterisationPolyconePhi::ComputeDimensions():" << G4endl;
347 pcone.DumpInfo();
348 }
349#endif
350}
351
352//---------------------------------------------------------------------
353G4ParameterisationPolyconeZ::
354G4ParameterisationPolyconeZ( EAxis axis, G4int nDiv,
355 G4double width, G4double offset,
356 G4VSolid* msolid, DivisionType divType)
357 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType )
358{
359
360 CheckParametersValidity();
361 SetType( "DivisionPolyconeZ" );
362
363 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
364 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
365
366 if( divType == DivWIDTH )
367 {
368 fnDiv =
369 CalculateNDiv( origparamMother->Z_values[origparamMother->Num_z_planes-1]
370 - origparamMother->Z_values[0] , width, offset );
371 }
372 else if( divType == DivNDIV )
373 {
374 fwidth =
375 CalculateNDiv( origparamMother->Z_values[origparamMother->Num_z_planes-1]
376 - origparamMother->Z_values[0] , nDiv, offset );
377 }
378
379#ifdef G4DIVDEBUG
380 if( verbose >= 1 )
381 {
382 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
383 << nDiv << G4endl
384 << " Offset " << foffset << " = " << offset << G4endl
385 << " Width " << fwidth << " = " << width << G4endl;
386 }
387#endif
388}
389
390//---------------------------------------------------------------------
391G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ()
392{
393}
394
395//------------------------------------------------------------------------
396G4double G4ParameterisationPolyconeZ::GetMaxParameter() const
397{
398 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
399 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
400 return std::abs (origparamMother->Z_values[origparamMother->Num_z_planes-1]
401 -origparamMother->Z_values[0]);
402}
403
404//---------------------------------------------------------------------
405void G4ParameterisationPolyconeZ::CheckParametersValidity()
406{
407 G4VDivisionParameterisation::CheckParametersValidity();
408
409 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
410
411 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH )
412 {
413 G4cerr << "WARNING - "
414 << "G4ParameterisationPolyconeZ::CheckParametersValidity()"
415 << G4endl
416 << " Solid " << msol->GetName() << G4endl
417 << " Division along Z will be done splitting in the "
418 << "defined z_planes." << G4endl
419 << " WIDTH will not be used !" << G4endl;
420 }
421
422 if( foffset != 0. )
423 {
424 G4cerr << "WARNING - "
425 << "G4ParameterisationPolyconeZ::CheckParametersValidity()"
426 << G4endl
427 << " Solid " << msol->GetName() << G4endl
428 << " Division along Z will be done splitting in the "
429 << "defined z_planes." << G4endl
430 << " OFFSET will not be used !" << G4endl;
431 }
432
433 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
434
435 if( origparamMother->Num_z_planes-1 != fnDiv )
436 {
437 G4cerr << "ERROR - "
438 << "G4ParameterisationPolyconeZ::CheckParametersValidity()"
439 << G4endl
440 << " Division along Z will be done splitting in the defined"
441 << G4endl
442 << " z_planes, i.e, the number of division would be :"
443 << " " << origparamMother->Num_z_planes-1
444 << " instead of " << fnDiv << " !"
445 << G4endl;
446 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
447 "IllegalConstruct", FatalException,
448 "Not supported configuration.");
449 }
450}
451
452//---------------------------------------------------------------------
453void
454G4ParameterisationPolyconeZ::
455ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
456{
457 G4Polycone* msol = (G4Polycone*)(GetMotherSolid());
458
459 //----- set translation: along Z axis
460 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
461 G4double posi = (origparamMother->Z_values[copyNo]
462 + origparamMother->Z_values[copyNo+1])/2;
463 G4ThreeVector origin(0.,0.,posi);
464 physVol->SetTranslation( origin );
465
466 //----- calculate rotation matrix: unit
467
468#ifdef G4DIVDEBUG
469 if( verbose >= 2 )
470 {
471 G4cout << " G4ParameterisationPolyconeZ - position: " << posi << G4endl
472 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
473 << " - fwidth: " << fwidth/deg << G4endl;
474 }
475#endif
476
477 ChangeRotMatrix( physVol );
478
479#ifdef G4DIVDEBUG
480 if( verbose >= 2 )
481 {
482 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeZ "
483 << copyNo << G4endl
484 << " Position: " << origin << " - Width: " << fwidth
485 << " - Axis: " << faxis << G4endl;
486 }
487#endif
488}
489
490//---------------------------------------------------------------------
491void
492G4ParameterisationPolyconeZ::
493ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
494 const G4VPhysicalVolume* ) const
495{
496 // only for mother number of planes = 2!!
497 //
498 G4Polycone* msol = (G4Polycone*)(fmotherSolid);
499
500 G4PolyconeHistorical* origparamMother = msol->GetOriginalParameters();
501 G4PolyconeHistorical origparam( *origparamMother );
502
503 G4double posi = (origparamMother->Z_values[copyNo]
504 + origparamMother->Z_values[copyNo+1])/2;
505
506 origparam.Num_z_planes = 2;
507 origparam.Z_values[0] = origparamMother->Z_values[copyNo] - posi;
508 origparam.Z_values[1] = origparamMother->Z_values[copyNo+1] - posi;
509 origparam.Rmin[0] = origparamMother->Rmin[copyNo];
510 origparam.Rmin[1] = origparamMother->Rmin[copyNo+1];
511 origparam.Rmax[0] = origparamMother->Rmax[copyNo];
512 origparam.Rmax[1] = origparamMother->Rmax[copyNo+1];
513
514 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
515 pcone.Reset(); // reset to new solid parameters
516
517#ifdef G4DIVDEBUG
518 if( verbose >= 2 )
519 {
520 G4cout << "G4ParameterisationPolyconeZ::ComputeDimensions()" << G4endl
521 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
522 pcone.DumpInfo();
523 }
524#endif
525}
Note: See TracBrowser for help on using the repository browser.