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

Last change on this file since 1242 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 22.9 KB
RevLine 
[831]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//
[1228]27// $Id: G4ParameterisationPolycone.cc,v 1.17 2009/05/18 19:30:29 ivana Exp $
28// GEANT4 tag $Name: geant4-09-03 $
[831]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)
[1228]357 : G4VParameterisationPolycone( axis, nDiv, width, offset, msolid, divType ),
358 fNSegment(0),
359 fOrigParamMother(((G4Polycone*)fmotherSolid)->GetOriginalParameters())
[831]360{
361
362 CheckParametersValidity();
363 SetType( "DivisionPolyconeZ" );
364
365 if( divType == DivWIDTH )
366 {
367 fnDiv =
[1228]368 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
369 - fOrigParamMother->Z_values[0] , width, offset );
[831]370 }
371 else if( divType == DivNDIV )
372 {
373 fwidth =
[1228]374 CalculateNDiv( fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
375 - fOrigParamMother->Z_values[0] , nDiv, offset );
[831]376 }
377
378#ifdef G4DIVDEBUG
379 if( verbose >= 1 )
380 {
381 G4cout << " G4ParameterisationPolyconeZ - # divisions " << fnDiv << " = "
382 << nDiv << G4endl
383 << " Offset " << foffset << " = " << offset << G4endl
384 << " Width " << fwidth << " = " << width << G4endl;
385 }
386#endif
387}
388
389//---------------------------------------------------------------------
390G4ParameterisationPolyconeZ::~G4ParameterisationPolyconeZ()
391{
392}
393
394//------------------------------------------------------------------------
[1228]395G4double G4ParameterisationPolyconeZ::GetR(G4double z,
396 G4double z1, G4double r1,
397 G4double z2, G4double r2) const
398{
399 // Linear parameterisation:
400 // r = az + b
401 // a = (r1 - r2)/(z1-z2)
402 // b = r1 - a*z1
403
404 return (r1-r2)/(z1-z2)*z + ( r1 - (r1-r2)/(z1-z2)*z1 ) ;
405}
406
407//------------------------------------------------------------------------
408G4double G4ParameterisationPolyconeZ::GetRmin(G4double z, G4int nseg) const
409{
410// Get Rmin in the given z position for the given polycone segment
411
412 return GetR(z,
413 fOrigParamMother->Z_values[nseg],
414 fOrigParamMother->Rmin[nseg],
415 fOrigParamMother->Z_values[nseg+1],
416 fOrigParamMother->Rmin[nseg+1]);
417}
418
419//------------------------------------------------------------------------
420G4double G4ParameterisationPolyconeZ::GetRmax(G4double z, G4int nseg) const
421{
422// Get Rmax in the given z position for the given polycone segment
423
424 return GetR(z,
425 fOrigParamMother->Z_values[nseg],
426 fOrigParamMother->Rmax[nseg],
427 fOrigParamMother->Z_values[nseg+1],
428 fOrigParamMother->Rmax[nseg+1]);
429}
430
431//------------------------------------------------------------------------
[831]432G4double G4ParameterisationPolyconeZ::GetMaxParameter() const
433{
[1228]434 return std::abs (fOrigParamMother->Z_values[fOrigParamMother->Num_z_planes-1]
435 -fOrigParamMother->Z_values[0]);
[831]436}
437
438//---------------------------------------------------------------------
439void G4ParameterisationPolyconeZ::CheckParametersValidity()
440{
441 G4VDivisionParameterisation::CheckParametersValidity();
[1228]442
443 // Division will be following the mother polycone segments
444 if( fDivisionType == DivNDIV ) {
445 if( fnDiv > fOrigParamMother->Num_z_planes-1 ) {
446 G4cerr << "ERROR - "
447 << "G4ParameterisationPolyconeZ::CheckParametersValidity()"
448 << G4endl
449 << " Division along Z will be done splitting in the defined"
450 << G4endl
451 << " z_planes, i.e, the number of division would be :"
452 << " " << fOrigParamMother->Num_z_planes-1
453 << " instead of " << fnDiv << " !"
454 << G4endl;
455 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
[831]456 "IllegalConstruct", FatalException,
457 "Not supported configuration.");
[1228]458 }
459 }
460
461 // Division will be done within one polycone segment
462 // with applying given width and offset
463 if( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
464 // Check if divided region does not span over more
465 // than one z segment
466
467 G4int isegstart = -1; // number of the segment containing start position
468 G4int isegend = -1; // number of the segment containing end position
469
470 if ( ! fReflectedSolid ) {
471 // The start/end position of the divided region
472 G4double zstart
473 = fOrigParamMother->Z_values[0] + foffset;
474 G4double zend
475 = fOrigParamMother->Z_values[0] + foffset + fnDiv* fwidth;
476
477 G4int counter = 0;
478 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
479 // first segment
480 if ( zstart >= fOrigParamMother->Z_values[counter] &&
481 zstart < fOrigParamMother->Z_values[counter+1] ) {
482 isegstart = counter;
483 }
484 // last segment
485 if ( zend > fOrigParamMother->Z_values[counter] &&
486 zend <= fOrigParamMother->Z_values[counter+1] ) {
487 isegend = counter;
488 }
489 ++counter;
490 }
491 }
492 else {
493 // The start/end position of the divided region
494 G4double zstart
495 = fOrigParamMother->Z_values[0] - foffset;
496 G4double zend
497 = fOrigParamMother->Z_values[0] - ( foffset + fnDiv* fwidth);
498
499 G4int counter = 0;
500 while ( isegend < 0 && counter < fOrigParamMother->Num_z_planes-1 ) {
501 // first segment
502 if ( zstart <= fOrigParamMother->Z_values[counter] &&
503 zstart > fOrigParamMother->Z_values[counter+1] ) {
504 isegstart = counter;
505 }
506 // last segment
507 if ( zend < fOrigParamMother->Z_values[counter] &&
508 zend >= fOrigParamMother->Z_values[counter+1] ) {
509 isegend = counter;
510 }
511 ++counter;
512 }
513 }
514
515
516 if ( isegstart != isegend ) {
517 G4cerr << "WARNING - "
518 << "G4ParameterisationPolyconeZ::CheckParametersValidity()"
519 << G4endl
520 << " Division with user defined width." << G4endl
521 << " Solid " << fmotherSolid->GetName() << G4endl
522 << " Divided region is not between two z planes."
523 << G4endl;
524
525 G4Exception("G4ParameterisationPolyconeZ::CheckParametersValidity()",
526 "IllegalConstruct", FatalException,
527 "Not supported configuration.");
528 }
529
530 fNSegment = isegstart;
531 }
[831]532}
533
534//---------------------------------------------------------------------
535void
536G4ParameterisationPolyconeZ::
537ComputeTransformation( const G4int copyNo, G4VPhysicalVolume* physVol) const
538{
[1228]539 if ( fDivisionType == DivNDIV ) {
540 // The position of the centre of copyNo-th mother polycone segment
541 G4double posi
542 = ( fOrigParamMother->Z_values[copyNo]
543 + fOrigParamMother->Z_values[copyNo+1])/2;
544 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
545 }
546
547 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
548 // The position of the centre of copyNo-th division
549 G4double posi = fOrigParamMother->Z_values[0];
550
551 if ( ! fReflectedSolid )
552 posi += foffset + (2*copyNo + 1) * fwidth/2.;
553 else
554 posi -= foffset + (2*copyNo + 1) * fwidth/2.;
555
556 physVol->SetTranslation( G4ThreeVector(0, 0, posi) );
557 }
[831]558
559 //----- calculate rotation matrix: unit
560
561#ifdef G4DIVDEBUG
562 if( verbose >= 2 )
563 {
564 G4cout << " G4ParameterisationPolyconeZ - position: " << posi << G4endl
565 << " copyNo: " << copyNo << " - foffset: " << foffset/deg
566 << " - fwidth: " << fwidth/deg << G4endl;
567 }
568#endif
569
570 ChangeRotMatrix( physVol );
571
572#ifdef G4DIVDEBUG
573 if( verbose >= 2 )
574 {
575 G4cout << std::setprecision(8) << " G4ParameterisationPolyconeZ "
576 << copyNo << G4endl
577 << " Position: " << origin << " - Width: " << fwidth
578 << " - Axis: " << faxis << G4endl;
579 }
580#endif
581}
582
583//---------------------------------------------------------------------
584void
585G4ParameterisationPolyconeZ::
586ComputeDimensions( G4Polycone& pcone, const G4int copyNo,
587 const G4VPhysicalVolume* ) const
588{
589
[1228]590 // Define division solid
591 G4PolyconeHistorical origparam;
592 G4int nz = 2;
593 origparam.Num_z_planes = nz;
594 origparam.Start_angle = fOrigParamMother->Start_angle;
595 origparam.Opening_angle = fOrigParamMother->Opening_angle;
[831]596
[1228]597 // Define division solid z sections
598 origparam.Z_values = new G4double[nz];
599 origparam.Rmin = new G4double[nz];
600 origparam.Rmax = new G4double[nz];
[831]601
[1228]602 if ( fDivisionType == DivNDIV ) {
603 // The position of the centre of copyNo-th mother polycone segment
604 G4double posi = (fOrigParamMother->Z_values[copyNo]
605 + fOrigParamMother->Z_values[copyNo+1])/2;
[831]606
[1228]607 origparam.Z_values[0] = fOrigParamMother->Z_values[copyNo] - posi;
608 origparam.Z_values[1] = fOrigParamMother->Z_values[copyNo+1] - posi;
609 origparam.Rmin[0] = fOrigParamMother->Rmin[copyNo];
610 origparam.Rmin[1] = fOrigParamMother->Rmin[copyNo+1];
611 origparam.Rmax[0] = fOrigParamMother->Rmax[copyNo];
612 origparam.Rmax[1] = fOrigParamMother->Rmax[copyNo+1];
613 }
614
615 if ( fDivisionType == DivNDIVandWIDTH || fDivisionType == DivWIDTH ) {
616 if ( ! fReflectedSolid ) {
617 origparam.Z_values[0] = - fwidth/2.;
618 origparam.Z_values[1] = fwidth/2.;
619
620 // The position of the centre of copyNo-th division
621 G4double posi
622 = fOrigParamMother->Z_values[0] + foffset + (2*copyNo + 1) * fwidth/2.;
623
624 // The first and last z sides z values
625 G4double zstart = posi - fwidth/2.;
626 G4double zend = posi + fwidth/2.;
627 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
628 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
629 origparam.Rmin[1] = GetRmin(zend, fNSegment);
630 origparam.Rmax[1] = GetRmax(zend, fNSegment);
631 }
632 else {
633 origparam.Z_values[0] = fwidth/2.;
634 origparam.Z_values[1] = - fwidth/2.;
635
636 // The position of the centre of copyNo-th division
637 G4double posi
638 = fOrigParamMother->Z_values[0] - ( foffset + (2*copyNo + 1) * fwidth/2.);
639
640 // The first and last z sides z values
641 G4double zstart = posi + fwidth/2.;
642 G4double zend = posi - fwidth/2.;
643 origparam.Rmin[0] = GetRmin(zstart, fNSegment);
644 origparam.Rmax[0] = GetRmax(zstart, fNSegment);
645 origparam.Rmin[1] = GetRmin(zend, fNSegment);
646 origparam.Rmax[1] = GetRmax(zend, fNSegment);
647 }
648
649 // It can happen due to rounding errors
650 if ( origparam.Rmin[0] < 0.0 ) origparam.Rmin[0] = 0.0;
651 if ( origparam.Rmin[nz-1] < 0.0 ) origparam.Rmin[1] = 0.0;
652 }
653
[831]654 pcone.SetOriginalParameters(&origparam); // copy values & transfer pointers
655 pcone.Reset(); // reset to new solid parameters
656
657#ifdef G4DIVDEBUG
658 if( verbose >= 2 )
659 {
660 G4cout << "G4ParameterisationPolyconeZ::ComputeDimensions()" << G4endl
661 << "-- Parametrised pcone copy-number: " << copyNo << G4endl;
662 pcone.DumpInfo();
663 }
664#endif
665}
Note: See TracBrowser for help on using the repository browser.