#include <BLCoordinates.hh>
Whenever an element is placed into the world group, the current centerline coordinates are normally used, rather than global coordinates (the local coordinates of the world volume/group).
Centerline coordinates are X,Y,Z, with Z along the nominal beam direction, X is beam left, and Y is beam up. The centerline itself is always X=0,Y=0, and consists of a sequence of segmentCLs each of which is straight, and the centerline is continuous at each corner (bend in the centerline). Except within bending magnets the reference particle is supposed to travel directly down the centerline. Z is continuous along the centerline, including at corners; the first segmentCL extends to Z=-infinity, and the last segmentCL extends to Z=+infinity.
Note that the coordinate transform global->centerline is NOT unique -- near corners both segmentCLs may refer to a given point. The transform centerline->global is unique for a given segmentCL. The centerline itself is continuous at corners, but tracks will in general have a discontinuity in their centerline coordinates when transitioning from one segmentCL to the next. A track will transition from one segmentCL to the next when either its Z coordinate exceeds that of the initial segmentCL's end (corner), or when its radius exceeds the radiusCut for the initial segmentCL.
Initially the centerline coordinates are identical to global coordinates. The start() function will change the starting location and orientation of the first segmentCL, and the corner() function ends the current segmentCL and begins a new one. start() must not be called after any element has been placed in the first segmentCL, or after corner() has been called.
Instances of this class represent the current centerline coordinates of a track, and are linked into the track via G4Track::SetUserInformation() (an instance is registered with BLManager to manage the instance linked into the track, and they can also be used for individual points, if necessary). Static members and functions are used to access and create the centerline segmentCLs.
Because of the need to transition tracks from one segmentCL to the next, corners with acute angles are prohibited.
NOTE: BLManager expects G4Track::GetUserInformation() to return a pointer to a BLCoordinates object, or to an object of a class derived from it. This means that any G4VUserTrackInformation class used in a simulation using BLManager should be derived from this class, and not directly from G4VUserTrackInformation (but verify is used to detect if it is not a BLCoordinates object).
Note the BLManager::UserSteppingAction() manages the updates to the current track's BLCoordinates before doing anything else.
NOTE: "CL" in a vairable means it is for centerline coordinates; "RC" means it is for reference coordinates.
Public Member Functions | |
BLCoordinates () | |
Default constructor -- for UserSteppingAction or without point. | |
BLCoordinates (const G4double _global[4]) | |
Constructor -- finds the earliest segmentCL containing the point (global coordinates). | |
virtual | ~BLCoordinates () |
Destructor. | |
void * | operator new (size_t) |
operator new - using G4Allocator. NOTE: any derived class MUST re-implement this. | |
void | operator delete (void *p) |
operator delete - using G4Allocator. NOTE: any derived class MUST re-implement this. | |
G4bool | isValid () const |
isValid() verifies that this is a valid BLCoordinates instance. NOTE: this is not virtual, and cannot be made virtual! | |
G4bool | setGlobal (const G4double _global[4]) |
setGlobal() updates the centerline coordinates for the new point[] of the track; uses the current segmentCL #, and will update it if the track crosses into a new segmentCL. returns false if the track fails radiusCut in its current and next segmentCL (or all segmentCL if ring is true). Also updates reference coordinates, if validRC is true and the centerline radiusCut succeeded. | |
G4bool | setGlobal (const G4ThreeVector &pos, G4double time) |
void | getCenterlineCoords (G4double local[4]) const |
getCenterlineCoords() returns the centerline coordinates. DEPRECATED: use getCoords() instead. | |
void | getCenterlineCoords (G4ThreeVector &local) const |
void | getCoords (BLCoordinateType type, G4double local[4]) |
getCoords() returns the appropriate coordinate values. | |
void | getCoords (BLCoordinateType type, G4ThreeVector &local) |
G4double | getCLZ () |
getCLZ() returns the centerline z coordinate. | |
G4RotationMatrix & | getRotation () |
getRotation() returns the appropriate coordinate rotation for the previous call to getCoords(). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: coord->getCoords(BLCOORD_CENTERLINE,local); B = coord->getRotation() * B; | |
G4RotationMatrix & | getRotation (BLCoordinateType type) |
getRotation() returns the rotation for the specified coordinates. Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = coord->getRotation(BLCOORD_CENTERLINE) * B; | |
G4RotationMatrix | getCenterlineRotation () |
getCenterlineRotation() returns the rotation matrix used for the last setGlobal() of this instance. This is the rotation of the centerline wrt global coordinate axes. DEPRECATED -- use getRotation(). | |
G4RotationMatrix | getTransformRotation () |
getTransformRotation() returns the rotation matrix of the global -> centerline coordinate transformation; it is the inverse of getCenterlineRotation(). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = getTransformRotation() * B; DEPRECATED -- use getRotation(). | |
int | getSegmentCL () |
getSegmentCL() re turns the CL segment # (for debugging). | |
virtual void | Print () const |
Print from G4VUserTrackInformation. | |
Static Public Member Functions | |
static double | getCurrentRadiusCut () |
getCurrentRadiusCut() returns the radiusCut for the current CL segment. | |
static void | update (G4Track *track) |
update() This function ensures that a UserTrackInformation object is linked into the track; if none is, it creates a new BLCoordinates instance and links it in. It also calls setGlobal() for the BLCoordinates linked into the track. | |
static void | start (const G4double point[4], G4double z, G4RotationMatrix *rotation=0, G4double radiusCut=0.0, G4bool _ring=false) |
start() sets the starting point and direction for the first segmentCL. Must not be called after corner() is called, or after any element has been placed using centerline coordinates. point[] is in global coordinates, z is the starting value of z in CENTERLINE coordinates, _ring indicates whether a ring is present (so jumps in segmentCL # are permitted), and rotation should be the rotation of the first segmentCL's centerline wrt the global coordinate axes (i.e. it is an active rotation of the centerline as an object wrt the global coordinate axes). | |
static void | corner (G4double z, G4RotationMatrix *rotation, G4double radiusCut) |
corner() creates a corner in the centerline by ending the current segmentCL and beginning a new one. z is the position in CENTERLINE coordinates of the current segmentCL, and rotation is relative to the current CENTERLINE coordinates (i.e. as an active rotation of the centerline as an object wrt the CENTERLINE coordinate axes of the current segmentCL). z must be greater than the centerline z of the previous corner or start. rotation can be an identity matrix (e.g. to merely set a new radiusCut). | |
static BLCoordinateTransform & | getCurrentTransform () |
getCurrentTransform() returns the current BLCoordinateTransform for the centerline coordinates (i.e. the coordinate transform global->centerline of the current segmentCL). | |
static G4RotationMatrix * | getCurrentRotation () |
getCurrentRotation() returns the current G4RotationMatrix for the centerline coordinates (i.e. the active rotation of the current segmentCL, as an object rotated wrt global coordinate axes). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = *getCurrentRotation() * B; | |
static void | getCurrentGlobal (const G4double local[4], G4double _global[4]) |
getCurrentGlobal() returns global coordinates from centerline coordinates, using the current segmentCL. | |
static void | getCurrentGlobal (const G4ThreeVector &local, G4ThreeVector &_global) |
getCurrentGlobal() returns global coordinates from centerline coordinates, using the current segmentCL. | |
static void | getGlobalAnywhere (const G4ThreeVector &local, G4ThreeVector &_global) |
getGlobalAnywhere() returns global coordinates from centerline coordinates, using the lowest segmentCL containing local. | |
static G4RotationMatrix * | getGlobalRotation (const G4ThreeVector &local) |
getGlobalRotation() returns the local->global rotation matrix coordinates, using the lowest segmentCL containing local. | |
static G4double | getRadiusCut (const G4ThreeVector &local) |
getRadiusCut() returns the radius cut for the lowest segmentCL containing local. | |
static BLCoordinateType | getCoordinateType (G4String s) |
getCoordinateType() converts a G4String name to a BLCoordinateType. Names are: global, centerline, reference; only the first letter is used, case insensitive. Also enables the generation of REference Coordinates, if they are selected by s. | |
static void | useReferenceCoordinates () |
useReferenceCoordinates() will arrange to fill segmentRCVector[] with reference coordinate segments as the reference particle is tracked. Must be called before the reference particle is tracked. Unless this is called, reference coordinates will not be available. Multiple calls are the same as one. | |
static void | setMostRecentReferenceZ (G4double z) |
setMostRecentReferenceZ() sets the most recent z value of the reference coords. VALID ONLY WHILE TRACKING REFERENCE PARTICLE! | |
static G4double | getMostRecentReferenceZ () |
getMostRecentReferenceZ() returns the most recent z value of the reference coords. VALID ONLY WHILE TRACKING REFERENCE PARTICLE! | |
Private Types | |
enum | { VERIFY = 0x13579246 } |
enum | State { UNUSED, CONSTRUCTING, TRACKING } |
Private Member Functions | |
G4bool | findCLSegment (const G4double _global[4]) |
void | updateReferenceCoordinates () |
updateReferenceCoordinates() will update the reference coordinates, given the centerline coordintes are already updated and valid. Called by setGlobal(). | |
Static Private Member Functions | |
static void | init () |
Private Attributes | |
int | verify |
unsigned int | segmentCL |
unsigned int | segmentRC |
BLCoordinateType | prevType |
G4double | global [4] |
G4double | centerline [4] |
G4double | reference [4] |
int | segmentCLRadiusFailure |
int | referenceCLsegment |
BLCoordinateTransform * | referenceTransformUsed |
Static Private Attributes | |
static G4Allocator< BLCoordinates > | allocator |
static std::vector< CLSegment > | segmentCLVector |
static unsigned | currentCL = 0 |
static G4bool | ring = false |
static State | state = UNUSED |
static G4bool | steppingVerbose = false |
static G4bool | useRC = false |
static G4bool | validRC = false |
static std::vector< RCSegment > | segmentRCVector |
static int | currentRC |
static G4double | tolerance = 1.0e-4 |
static G4double | mostRecentReferenceZ = 0 |
Classes | |
struct | CLSegment |
struct | RCSegment |
class | ReferenceCoordinates |
anonymous enum [private] |
enum BLCoordinates::State [private] |
BLCoordinates::BLCoordinates | ( | ) |
Default constructor -- for UserSteppingAction or without point.
References centerline, global, reference, referenceCLsegment, referenceTransformUsed, segmentCL, segmentCLRadiusFailure, segmentRC, VERIFY, and verify.
Referenced by update().
00065 : G4VUserTrackInformation() 00066 { 00067 // no state change, yet 00068 verify=VERIFY; 00069 segmentCL=99999999; 00070 segmentRC=99999999; 00071 referenceCLsegment = -1; 00072 referenceTransformUsed = 0; 00073 global[0]=global[1]=global[2]=global[3]=0.0; 00074 centerline[0]=centerline[1]=centerline[2]=centerline[3]=0.0; 00075 reference[0]=reference[1]=reference[2]=reference[3]=0.0; 00076 segmentCLRadiusFailure = -1; 00077 }
BLCoordinates::BLCoordinates | ( | const G4double | _global[4] | ) |
Constructor -- finds the earliest segmentCL containing the point (global coordinates).
References findCLSegment(), state, TRACKING, VERIFY, and verify.
00079 : 00080 G4VUserTrackInformation() 00081 { 00082 state = TRACKING; 00083 verify = VERIFY; 00084 findCLSegment(_global); 00085 }
G4bool BLCoordinates::findCLSegment | ( | const G4double | _global[4] | ) | [private] |
References centerline, segmentCL, segmentCLRadiusFailure, and segmentCLVector.
Referenced by BLCoordinates(), and setGlobal().
00094 { 00095 for(segmentCL=0; segmentCL<segmentCLVector.size(); ++segmentCL) { 00096 segmentCLVector[segmentCL].transform.getLocal(centerline,_global); 00097 //if(steppingVerbose) printf("FindCLSegment: seg=%d minZ=%.3f maxZ=%.3f Zcl=%.3f Rcl=%.3f rCut=%.3f\n",segmentCL,segmentCLVector[segmentCL].minZ,segmentCLVector[segmentCL].maxZ,centerline[2],sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]),segmentCLVector[segmentCL].radiusCut); 00098 // note no check against minZ -- segmentCLs are ordered 00099 if(centerline[2] <= segmentCLVector[segmentCL].maxZ && 00100 (segmentCLVector[segmentCL].radiusCut <= 0.0 || 00101 (int)segmentCL == segmentCLRadiusFailure || 00102 sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]) <= 00103 segmentCLVector[segmentCL].radiusCut)) { 00104 return true; 00105 } 00106 } 00107 --segmentCL; 00108 return false; 00109 }
void* BLCoordinates::operator new | ( | size_t | ) | [inline] |
operator new - using G4Allocator. NOTE: any derived class MUST re-implement this.
Reimplemented in BLTrackInfo.
References allocator.
00119 { 00120 return (void *)allocator.MallocSingle(); 00121 }
void BLCoordinates::operator delete | ( | void * | p | ) | [inline] |
operator delete - using G4Allocator. NOTE: any derived class MUST re-implement this.
Reimplemented in BLTrackInfo.
References allocator.
00125 { 00126 allocator.FreeSingle((BLCoordinates *)p); 00127 }
G4bool BLCoordinates::isValid | ( | ) | const |
isValid() verifies that this is a valid BLCoordinates instance. NOTE: this is not virtual, and cannot be made virtual!
References VERIFY, and verify.
Referenced by BLTrackNTuple::appendTrack(), BLCMDmovie::beamStep(), BLManager::getExternalParentID(), BLManager::getExternalTrackID(), BLManager::PreUserTrackingAction(), BLCoordinates::ReferenceCoordinates::PreUserTrackingAction(), BLCMDmovie::referenceStep(), BLManager::setExternalTrackID(), BLEvaluator::setTrackVariables(), BLManager::steppingVerbosePrint(), update(), BLManager::UserSteppingAction(), TrackerPlaneInstance::UserSteppingAction(), BLCMDtracker::UserZSteppingAction(), and BLCMDprofile::Entry::UserZSteppingAction().
G4bool BLCoordinates::setGlobal | ( | const G4double | _global[4] | ) |
setGlobal() updates the centerline coordinates for the new point[] of the track; uses the current segmentCL #, and will update it if the track crosses into a new segmentCL. returns false if the track fails radiusCut in its current and next segmentCL (or all segmentCL if ring is true). Also updates reference coordinates, if validRC is true and the centerline radiusCut succeeded.
References centerline, findCLSegment(), global, ring, segmentCL, segmentCLRadiusFailure, segmentCLVector, updateReferenceCoordinates(), and validRC.
Referenced by BLCMDmovie::dumpGeometry(), BLManager::PreUserTrackingAction(), setGlobal(), and update().
00112 { 00113 global[0] = _global[0]; 00114 global[1] = _global[1]; 00115 global[2] = _global[2]; 00116 global[3] = _global[3]; 00117 00118 if(segmentCL >= segmentCLVector.size()) { 00119 //printf("setGlobal-0: calling findCLSegment\n"); 00120 G4bool ret = findCLSegment(_global); 00121 if(ret && validRC) updateReferenceCoordinates(); 00122 //if(!ret) printf("setGlobal-0: findCLSegment FAILED\n"); 00123 return ret; 00124 } 00125 00126 segmentCLVector[segmentCL].transform.getLocal(centerline,_global); 00127 //printf("setGlobal-1A: global=%.1f,%.1f,%.1f cl=%.1f,%.1f,%.1f\n",_global[0],_global[1],_global[2],centerline[0],centerline[1],centerline[2]); 00128 //printf("setGlobal-1: seg=%d minZ=%.3f maxZ=%.3f Zcl=%.3f Rcl=%.3f rCut=%.3f\n",segmentCL,segmentCLVector[segmentCL].minZ,segmentCLVector[segmentCL].maxZ,centerline[2],sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]),segmentCLVector[segmentCL].radiusCut); 00129 while(centerline[2] > segmentCLVector[segmentCL].maxZ || 00130 (segmentCLVector[segmentCL].radiusCut > 0.0 && 00131 (int)segmentCL != segmentCLRadiusFailure && 00132 sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]) > 00133 segmentCLVector[segmentCL].radiusCut)) { 00134 if(++segmentCL >= segmentCLVector.size()) { 00135 if(ring && findCLSegment(_global)) { 00136 segmentCLVector[segmentCL].transform.getLocal(centerline, 00137 _global); 00138 if(validRC) updateReferenceCoordinates(); 00139 return true; 00140 } 00141 --segmentCL; 00142 //printf("setGlobal-1a: FAILED\n"); 00143 return false; 00144 } 00145 segmentCLVector[segmentCL].transform.getLocal(centerline,_global); 00146 //printf("setGlobal-2: seg=%d minZ=%.3f maxZ=%.3f Zcl=%.3f Rcl=%.3f rCut=%.3f\n",segmentCL,segmentCLVector[segmentCL].minZ,segmentCLVector[segmentCL].maxZ,centerline[2],sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]),segmentCLVector[segmentCL].radiusCut); 00147 /*** 00148 // note no check against minZ -- segmentCLs are ordered 00149 if(centerline[2] > segmentCLVector[segmentCL].maxZ || 00150 (segmentCLVector[segmentCL].radiusCut > 0.0 && 00151 segmentCL != segmentCLRadiusFailure && 00152 sqrt(centerline[0]*centerline[0]+centerline[1]*centerline[1]) > 00153 segmentCLVector[segmentCL].radiusCut)) { 00154 if(ring && findCLSegment(_global)) { 00155 segmentCLVector[segmentCL].transform.getLocal(centerline, 00156 _global); 00157 if(validRC) updateReferenceCoordinates(); 00158 return true; 00159 } 00160 //printf("setGlobal-2a: FAILED\n"); 00161 return false; 00162 } 00163 ***/ 00164 } 00165 00166 if(validRC) updateReferenceCoordinates(); 00167 return true; 00168 }
G4bool BLCoordinates::setGlobal | ( | const G4ThreeVector & | pos, | |
G4double | time | |||
) |
References global, and setGlobal().
00171 { 00172 global[0]=pos[0],global[1]=pos[1],global[2]=pos[2]; 00173 global[3]=time; 00174 return setGlobal(global); 00175 }
void BLCoordinates::getCenterlineCoords | ( | G4double | local[4] | ) | const [inline] |
getCenterlineCoords() returns the centerline coordinates. DEPRECATED: use getCoords() instead.
References centerline.
00146 { local[0]=centerline[0]; local[1]=centerline[1]; 00147 local[2]=centerline[2]; local[3]=centerline[3]; }
void BLCoordinates::getCenterlineCoords | ( | G4ThreeVector & | local | ) | const [inline] |
References centerline.
00149 { local[0]=centerline[0]; local[1]=centerline[1]; 00150 local[2]=centerline[2]; }
void BLCoordinates::getCoords | ( | BLCoordinateType | type, | |
G4double | local[4] | |||
) |
getCoords() returns the appropriate coordinate values.
References BLAssert, BLCOORD_CENTERLINE, BLCOORD_GLOBAL, BLCOORD_REFERENCE, centerline, BLCoordinateTransform::getLocal(), BLCoordinateTransform::getRotation(), global, prevType, reference, referenceTransformUsed, segmentCL, segmentCLVector, segmentRC, segmentRCVector, and validRC.
Referenced by BLTrackNTuple::appendTrack(), BLCMDmovie::beamStep(), BLCMDmovie::dumpGeometry(), getCoords(), BLCoordinates::ReferenceCoordinates::PreUserTrackingAction(), BLCMDmovie::referenceStep(), BLEvaluator::setTrackVariables(), BLManager::steppingVerbosePrint(), TrackerPlaneInstance::UserSteppingAction(), BLCMDtune::UserZSteppingAction(), BLCMDtracker::UserZSteppingAction(), and BLCMDprofile::Entry::UserZSteppingAction().
00196 { 00197 prevType = type; 00198 00199 switch(type) { 00200 case BLCOORD_GLOBAL: 00201 local[0] = global[0]; local[1] = global[1]; 00202 local[2] = global[2]; local[3] = global[3]; 00203 //printf("getCoords(GLOBAL): %.1f,%.1f,%.1f\n",local[0],local[1],local[2]); 00204 return; 00205 case BLCOORD_CENTERLINE: 00206 BLAssert(segmentCL < segmentCLVector.size()); 00207 local[0] = centerline[0]; local[1] = centerline[1]; 00208 local[2] = centerline[2]; local[3] = centerline[3]; 00209 //printf("getCoords(CL): %.1f,%.1f,%.1f\n",local[0],local[1],local[2]); 00210 return; 00211 case BLCOORD_REFERENCE: 00212 BLAssert(validRC); 00213 // interpolate 00214 if(reference[2] < segmentRCVector[segmentRC].interpolatePrev && 00215 segmentRC > 0) { 00216 referenceTransformUsed = 00217 &segmentRCVector[segmentRC-1].interpolate; 00218 referenceTransformUsed->getLocal(reference,global); 00219 } else if(reference[2] > 00220 segmentRCVector[segmentRC].interpolateThis) { 00221 referenceTransformUsed = 00222 &segmentRCVector[segmentRC].interpolate; 00223 referenceTransformUsed->getLocal(reference,global); 00224 } 00225 G4ThreeVector axis = referenceTransformUsed->getRotation().getAxis(); 00226 G4double angle = referenceTransformUsed->getRotation().getDelta(); 00227 //printf("segmentRC=%d minZ=%.4g maxZ=%.4g global=%.1f,%.1f,%.1f ref=%.1f,%.1f,%.1f axis=%.1f,%.1f,%.1f angle(deg)=%.1f\n",segmentRC,segmentRCVector[segmentRC].minZ,segmentRCVector[segmentRC].maxZ,global[0],global[1],global[2],reference[0],reference[1],reference[2],axis[0],axis[1],axis[2],angle/deg); 00228 local[0] = reference[0]; local[1] = reference[1]; 00229 local[2] = reference[2]; local[3] = reference[3]; 00230 //printf("getCoords(REF): %.1f,%.1f,%.1f\n",local[0],local[1],local[2]); 00231 return; 00232 } 00233 G4Exception("BLCoordinates","Invalid Coordinate Type",FatalException, 00234 ""); 00235 }
void BLCoordinates::getCoords | ( | BLCoordinateType | type, | |
G4ThreeVector & | local | |||
) |
References getCoords().
00238 { 00239 G4double tmp[4]; 00240 getCoords(type,tmp); 00241 local[0] = tmp[0]; 00242 local[1] = tmp[1]; 00243 local[2] = tmp[2]; 00244 //printf("getCoords(): %.1f,%.1f,%.1f\n",local[0],local[1],local[2]); 00245 }
G4double BLCoordinates::getCLZ | ( | ) | [inline] |
getCLZ() returns the centerline z coordinate.
References centerline.
Referenced by BLManager::PreUserTrackingAction(), and BLManager::UserSteppingAction().
00157 { return centerline[2]; }
G4RotationMatrix& BLCoordinates::getRotation | ( | ) | [inline] |
getRotation() returns the appropriate coordinate rotation for the previous call to getCoords(). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: coord->getCoords(BLCOORD_CENTERLINE,local); B = coord->getRotation() * B;
References prevType.
Referenced by BLTrackNTuple::appendTrack(), BLCMDmovie::beamStep(), BLEvaluator::setTrackVariables(), BLManager::steppingVerbosePrint(), TrackerPlaneInstance::UserSteppingAction(), BLCMDtracker::UserZSteppingAction(), and BLCMDprofile::Entry::UserZSteppingAction().
00165 { return getRotation(prevType); }
G4RotationMatrix & BLCoordinates::getRotation | ( | BLCoordinateType | type | ) |
getRotation() returns the rotation for the specified coordinates. Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = coord->getRotation(BLCOORD_CENTERLINE) * B;
References BLAssert, BLCOORD_CENTERLINE, BLCOORD_GLOBAL, BLCOORD_REFERENCE, BLCoordinateTransform::getRotation(), IdentityMatrix, referenceTransformUsed, segmentCL, segmentCLVector, and validRC.
00248 { 00249 switch(type) { 00250 case BLCOORD_GLOBAL: 00251 return IdentityMatrix; 00252 case BLCOORD_CENTERLINE: 00253 BLAssert(segmentCL < segmentCLVector.size()); 00254 return segmentCLVector[segmentCL].transform.getRotation(); 00255 case BLCOORD_REFERENCE: 00256 BLAssert(validRC && referenceTransformUsed != 0); 00257 return referenceTransformUsed->getRotation(); 00258 } 00259 G4Exception("BLCoordinates","Invalid Coordinate Type",FatalException, 00260 ""); 00261 return IdentityMatrix; 00262 }
G4RotationMatrix BLCoordinates::getCenterlineRotation | ( | ) | [inline] |
getCenterlineRotation() returns the rotation matrix used for the last setGlobal() of this instance. This is the rotation of the centerline wrt global coordinate axes. DEPRECATED -- use getRotation().
References segmentCL, and segmentCLVector.
00177 { 00178 if(segmentCL >= segmentCLVector.size()) 00179 return G4RotationMatrix(); 00180 return segmentCLVector[segmentCL].transform.getRotationInverse(); 00181 }
static double BLCoordinates::getCurrentRadiusCut | ( | ) | [inline, static] |
getCurrentRadiusCut() returns the radiusCut for the current CL segment.
References currentCL, and segmentCLVector.
Referenced by BLCMDcorner::command().
00185 { 00186 if(currentCL >= segmentCLVector.size()) 00187 return DBL_MAX; 00188 else 00189 return segmentCLVector[currentCL].radiusCut; 00190 }
G4RotationMatrix BLCoordinates::getTransformRotation | ( | ) | [inline] |
getTransformRotation() returns the rotation matrix of the global -> centerline coordinate transformation; it is the inverse of getCenterlineRotation(). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = getTransformRotation() * B; DEPRECATED -- use getRotation().
References segmentCL, and segmentCLVector.
00199 { 00200 if(segmentCL >= segmentCLVector.size()) 00201 return G4RotationMatrix(); 00202 return segmentCLVector[segmentCL].transform.getRotation(); 00203 }
int BLCoordinates::getSegmentCL | ( | ) | [inline] |
getSegmentCL() re turns the CL segment # (for debugging).
References segmentCL.
Referenced by BLManager::steppingVerbosePrint().
00206 { return segmentCL; }
void BLCoordinates::update | ( | G4Track * | track | ) | [static] |
update() This function ensures that a UserTrackInformation object is linked into the track; if none is, it creates a new BLCoordinates instance and links it in. It also calls setGlobal() for the BLCoordinates linked into the track.
References BLCoordinates(), isValid(), segmentCL, segmentCLRadiusFailure, setGlobal(), and steppingVerbose.
Referenced by BLManager::UserSteppingAction().
00178 { 00179 // ensure there is a BLCoordiantes linked into the track 00180 BLCoordinates *c = (BLCoordinates *)track->GetUserInformation(); 00181 if(!c || !c->isValid()) { 00182 c = new BLCoordinates(); 00183 track->SetUserInformation(c); 00184 } 00185 00186 // now update the object and apply radiusCut 00187 if(!c->setGlobal(track->GetPosition(),track->GetGlobalTime())) { 00188 track->SetTrackStatus(fStopAndKill); 00189 if(steppingVerbose) 00190 printf("Track Centerline radius exceeds radiusCut\n"); 00191 c->segmentCLRadiusFailure = c->segmentCL; 00192 } 00193 }
virtual void BLCoordinates::Print | ( | ) | const [inline, virtual] |
void BLCoordinates::init | ( | ) | [static, private] |
References CONSTRUCTING, currentCL, BLManager::getObject(), BLManager::getSteppingVerbose(), segmentCLVector, state, steppingVerbose, TRACKING, and UNUSED.
Referenced by corner(), getCurrentGlobal(), getCurrentRotation(), getCurrentTransform(), and start().
00265 { 00266 if(state == TRACKING) { 00267 G4Exception("BLCoordinates","Invalid init during tracking", 00268 FatalException, ""); 00269 } 00270 if(state == UNUSED) { 00271 state = CONSTRUCTING; 00272 segmentCLVector.push_back(CLSegment(-DBL_MAX, +DBL_MAX, 0.0, 00273 BLCoordinateTransform())); 00274 currentCL = 0; 00275 steppingVerbose = 00276 BLManager::getObject()->getSteppingVerbose() != 0; 00277 } 00278 }
void BLCoordinates::start | ( | const G4double | point[4], | |
G4double | z, | |||
G4RotationMatrix * | rotation = 0 , |
|||
G4double | radiusCut = 0.0 , |
|||
G4bool | _ring = false | |||
) | [static] |
start() sets the starting point and direction for the first segmentCL. Must not be called after corner() is called, or after any element has been placed using centerline coordinates. point[] is in global coordinates, z is the starting value of z in CENTERLINE coordinates, _ring indicates whether a ring is present (so jumps in segmentCL # are permitted), and rotation should be the rotation of the first segmentCL's centerline wrt the global coordinate axes (i.e. it is an active rotation of the centerline as an object wrt the global coordinate axes).
References currentCL, BLManager::getObject(), BLManager::getSteppingVerbose(), init(), ring, segmentCLVector, state, steppingVerbose, and UNUSED.
Referenced by BLCMDstart::command().
00283 { 00284 if(state != UNUSED || currentCL != 0) { 00285 G4Exception("BLCoordinates","Invalid start",FatalException, ""); 00286 } 00287 init(); 00288 00289 G4ThreeVector pos(0,0,-z); 00290 if(rotation) 00291 pos = *rotation * pos; 00292 pos[0] += _global[0]; 00293 pos[1] += _global[1]; 00294 pos[2] += _global[2]; 00295 segmentCLVector[0].transform = (rotation ? 00296 BLCoordinateTransform(rotation,pos) : 00297 BLCoordinateTransform(pos)); 00298 segmentCLVector[0].radiusCut = radiusCut; 00299 ring = _ring; 00300 steppingVerbose = BLManager::getObject()->getSteppingVerbose() != 0; 00301 }
void BLCoordinates::corner | ( | G4double | z, | |
G4RotationMatrix * | rotation, | |||
G4double | radiusCut | |||
) | [static] |
corner() creates a corner in the centerline by ending the current segmentCL and beginning a new one. z is the position in CENTERLINE coordinates of the current segmentCL, and rotation is relative to the current CENTERLINE coordinates (i.e. as an active rotation of the centerline as an object wrt the CENTERLINE coordinate axes of the current segmentCL). z must be greater than the centerline z of the previous corner or start. rotation can be an identity matrix (e.g. to merely set a new radiusCut).
References currentCL, getCurrentRotation(), getCurrentTransform(), init(), and segmentCLVector.
Referenced by BLCMDcornerarc::command(), and BLCMDcorner::command().
00305 { 00306 init(); 00307 00308 G4double local[4], glob[4]; 00309 local[0]=local[1]=local[3]=0.0; 00310 local[2] = z; 00311 segmentCLVector[currentCL].transform.getGlobal(local,glob); 00312 00313 segmentCLVector[currentCL].maxZ = z; 00314 G4RotationMatrix *p=0; 00315 if(rotation) 00316 // if C=currentCL rot, and R=this rot, we want R referenced to 00317 // the CURRENT centerline coordinates, so we need C R C^-1. 00318 // But we compose this with C (on the right), and get C R. 00319 p = new G4RotationMatrix(*getCurrentRotation() * *rotation); 00320 else if(getCurrentTransform().isRotated()) 00321 p = getCurrentRotation(); 00322 00323 G4ThreeVector pos(0,0,-z); 00324 if(p) 00325 pos = *p * pos; 00326 pos[0] += glob[0]; 00327 pos[1] += glob[1]; 00328 pos[2] += glob[2]; 00329 segmentCLVector.push_back(CLSegment(z, +DBL_MAX, radiusCut, 00330 BLCoordinateTransform(p,pos))); 00331 currentCL = segmentCLVector.size() - 1; 00332 }
BLCoordinateTransform & BLCoordinates::getCurrentTransform | ( | ) | [static] |
getCurrentTransform() returns the current BLCoordinateTransform for the centerline coordinates (i.e. the coordinate transform global->centerline of the current segmentCL).
References currentCL, init(), and segmentCLVector.
Referenced by corner().
00335 { 00336 init(); 00337 return segmentCLVector[currentCL].transform; 00338 }
G4RotationMatrix * BLCoordinates::getCurrentRotation | ( | ) | [static] |
getCurrentRotation() returns the current G4RotationMatrix for the centerline coordinates (i.e. the active rotation of the current segmentCL, as an object rotated wrt global coordinate axes). Given a G4ThreeVector B in global coordinates, this transforms it to centerline components: B = *getCurrentRotation() * B;
References currentCL, init(), and segmentCLVector.
Referenced by BLCMDtracker::command(), BLCMDreference::command(), BLCMDcorner::command(), BLCMDbeam::command(), corner(), and BLGroup::placeElement().
00341 { 00342 init(); 00343 return &segmentCLVector[currentCL].transform.getRotationInverse(); 00344 }
void BLCoordinates::getCurrentGlobal | ( | const G4double | local[4], | |
G4double | _global[4] | |||
) | [static] |
getCurrentGlobal() returns global coordinates from centerline coordinates, using the current segmentCL.
References currentCL, init(), and segmentCLVector.
Referenced by BLCMDtracker::command(), BLCMDreference::command(), BLCMDbeam::command(), and BLGroup::placeElement().
00347 { 00348 init(); 00349 segmentCLVector[currentCL].transform.getGlobal(local,_global); 00350 }
void BLCoordinates::getCurrentGlobal | ( | const G4ThreeVector & | local, | |
G4ThreeVector & | _global | |||
) | [static] |
getCurrentGlobal() returns global coordinates from centerline coordinates, using the current segmentCL.
References currentCL, init(), and segmentCLVector.
00353 { 00354 init(); 00355 segmentCLVector[currentCL].transform.getGlobal(local,_global); 00356 }
void BLCoordinates::getGlobalAnywhere | ( | const G4ThreeVector & | local, | |
G4ThreeVector & | _global | |||
) | [static] |
getGlobalAnywhere() returns global coordinates from centerline coordinates, using the lowest segmentCL containing local.
References segmentCLVector.
00360 { 00361 for(unsigned int seg=0; seg<segmentCLVector.size(); ++seg) { 00362 // note no check against minZ -- segmentCLs are ordered 00363 if(local[2] <= segmentCLVector[seg].maxZ && 00364 (segmentCLVector[seg].radiusCut <= 0.0 || 00365 sqrt(local[0]*local[0]+local[1]*local[1]) <= 00366 segmentCLVector[seg].radiusCut)) { 00367 segmentCLVector[seg].transform.getGlobal(local,_global); 00368 return; 00369 } 00370 } 00371 }
G4RotationMatrix * BLCoordinates::getGlobalRotation | ( | const G4ThreeVector & | local | ) | [static] |
getGlobalRotation() returns the local->global rotation matrix coordinates, using the lowest segmentCL containing local.
References segmentCLVector.
00374 { 00375 for(unsigned int seg=0; seg<segmentCLVector.size(); ++seg) { 00376 // note no check against minZ -- segmentCLs are ordered 00377 if(local[2] <= segmentCLVector[seg].maxZ && 00378 (segmentCLVector[seg].radiusCut <= 0.0 || 00379 sqrt(local[0]*local[0]+local[1]*local[1]) <= 00380 segmentCLVector[seg].radiusCut)) { 00381 return &segmentCLVector[seg].transform.getRotation(); 00382 } 00383 } 00384 00385 return 0; 00386 }
G4double BLCoordinates::getRadiusCut | ( | const G4ThreeVector & | local | ) | [static] |
getRadiusCut() returns the radius cut for the lowest segmentCL containing local.
References segmentCLVector.
00389 { 00390 for(unsigned int seg=0; seg<segmentCLVector.size(); ++seg) { 00391 // note no check against minZ -- segmentCLs are ordered 00392 if(local[2] <= segmentCLVector[seg].maxZ && 00393 (segmentCLVector[seg].radiusCut <= 0.0 || 00394 sqrt(local[0]*local[0]+local[1]*local[1]) <= 00395 segmentCLVector[seg].radiusCut)) { 00396 return segmentCLVector[seg].radiusCut; 00397 } 00398 } 00399 00400 return 0; 00401 }
BLCoordinateType BLCoordinates::getCoordinateType | ( | G4String | s | ) | [static] |
getCoordinateType() converts a G4String name to a BLCoordinateType. Names are: global, centerline, reference; only the first letter is used, case insensitive. Also enables the generation of REference Coordinates, if they are selected by s.
References BLCOORD_CENTERLINE, BLCOORD_GLOBAL, BLCOORD_REFERENCE, BLCommand::printError(), and useReferenceCoordinates().
Referenced by BLCMDzntuple::command(), BLCMDvirtualdetector::command(), BLCMDtrace::command(), BLCMDtimentuple::command(), BLCMDprofile::command(), BLCMDprintf::command(), BLCMDplace::command(), BLCMDnewparticlentuple::command(), BLCMDmovie::command(), and BLCMDbeamlossntuple::command().
00404 { 00405 char c = tolower(s[0u]); 00406 switch(c) { 00407 case 'g': return BLCOORD_GLOBAL; 00408 case 'c': return BLCOORD_CENTERLINE; 00409 case 'r': useReferenceCoordinates(); 00410 return BLCOORD_REFERENCE; 00411 default: BLCommand::printError("Unknown coordinate type '%s'," 00412 " Centerline used.\n",s.c_str()); 00413 return BLCOORD_CENTERLINE; 00414 } 00415 }
void BLCoordinates::updateReferenceCoordinates | ( | ) | [private] |
updateReferenceCoordinates() will update the reference coordinates, given the centerline coordintes are already updated and valid. Called by setGlobal().
References global, reference, referenceTransformUsed, segmentCL, segmentCLVector, segmentRC, and segmentRCVector.
Referenced by setGlobal().
00428 { 00429 // if we are lost, use CL coord segment as a starting point 00430 if(segmentRC >= segmentRCVector.size()) { 00431 segmentRC = segmentCLVector[segmentCL].firstRCSegment; 00432 if(segmentRC > 0) --segmentRC; // try 1 previous 00433 } 00434 if(segmentRC >= segmentRCVector.size()) 00435 segmentRC = 0; 00436 00437 for( ; segmentRC<segmentRCVector.size(); ++segmentRC) { 00438 segmentRCVector[segmentRC].transform.getLocal(reference,global); 00439 //printf("segmentRC=%d minZ=%.4g maxZ=%.4g global=%.1f,%.1f,%.1f ref=%.1f,%.1f,%.3f\n",segmentRC,segmentRCVector[segmentRC].minZ,segmentRCVector[segmentRC].maxZ,global[0],global[1],global[2],reference[0],reference[1],reference[2]); 00440 referenceTransformUsed = &segmentRCVector[segmentRC].transform; 00441 if(segmentRC != 0) // no correction (minZ = -DBL_MAX) 00442 reference[2] += (segmentRCVector[segmentRC].correction - 1.0) 00443 * (reference[2]-segmentRCVector[segmentRC].minZ); 00444 //printf(" corrected reference[2]=%.1f\n",reference[2]); 00445 if(reference[2] > segmentRCVector[segmentRC].maxZ) 00446 continue; 00447 return; 00448 } 00449 G4Exception("BLCoordinates","Cannot determine reference coordinate segment", 00450 FatalException,""); 00451 }
void BLCoordinates::useReferenceCoordinates | ( | ) | [static] |
useReferenceCoordinates() will arrange to fill segmentRCVector[] with reference coordinate segments as the reference particle is tracked. Must be called before the reference particle is tracked. Unless this is called, reference coordinates will not be available. Multiple calls are the same as one.
References useRC.
Referenced by BLCMDmovie::command(), and getCoordinateType().
00418 { 00419 00420 if(!useRC) { 00421 useRC = true; 00422 new ReferenceCoordinates(); 00423 } 00424 }
static void BLCoordinates::setMostRecentReferenceZ | ( | G4double | z | ) | [inline, static] |
setMostRecentReferenceZ() sets the most recent z value of the reference coords. VALID ONLY WHILE TRACKING REFERENCE PARTICLE!
References mostRecentReferenceZ.
Referenced by BLCoordinates::ReferenceCoordinates::PreUserTrackingAction(), and BLCoordinates::ReferenceCoordinates::UserSteppingAction().
00385 { mostRecentReferenceZ = z; }
static G4double BLCoordinates::getMostRecentReferenceZ | ( | ) | [inline, static] |
getMostRecentReferenceZ() returns the most recent z value of the reference coords. VALID ONLY WHILE TRACKING REFERENCE PARTICLE!
References mostRecentReferenceZ.
Referenced by BLCMDmovie::referenceStep().
00390 { return mostRecentReferenceZ; }
int BLCoordinates::verify [private] |
Referenced by BLCoordinates(), and isValid().
unsigned int BLCoordinates::segmentCL [private] |
unsigned int BLCoordinates::segmentRC [private] |
BLCoordinateType BLCoordinates::prevType [private] |
Referenced by getCoords(), and getRotation().
G4double BLCoordinates::global[4] [private] |
G4double BLCoordinates::centerline[4] [private] |
Referenced by BLCoordinates(), findCLSegment(), getCenterlineCoords(), getCLZ(), getCoords(), and setGlobal().
G4double BLCoordinates::reference[4] [private] |
Referenced by BLCoordinates(), getCoords(), and updateReferenceCoordinates().
int BLCoordinates::segmentCLRadiusFailure [private] |
Referenced by BLCoordinates(), findCLSegment(), setGlobal(), and update().
G4Allocator< BLCoordinates > BLCoordinates::allocator [static, private] |
Referenced by operator delete(), and operator new().
std::vector< BLCoordinates::CLSegment > BLCoordinates::segmentCLVector [static, private] |
Referenced by BLCoordinates::ReferenceCoordinates::addStep(), corner(), findCLSegment(), getCenterlineRotation(), getCoords(), getCurrentGlobal(), getCurrentRadiusCut(), getCurrentRotation(), getCurrentTransform(), getGlobalAnywhere(), getGlobalRotation(), getRadiusCut(), getRotation(), getTransformRotation(), init(), setGlobal(), start(), and updateReferenceCoordinates().
unsigned BLCoordinates::currentCL = 0 [static, private] |
G4bool BLCoordinates::ring = false [static, private] |
Referenced by setGlobal(), and start().
BLCoordinates::State BLCoordinates::state = UNUSED [static, private] |
Referenced by BLCoordinates(), init(), and start().
G4bool BLCoordinates::steppingVerbose = false [static, private] |
G4bool BLCoordinates::useRC = false [static, private] |
Referenced by useReferenceCoordinates().
G4bool BLCoordinates::validRC = false [static, private] |
std::vector< BLCoordinates::RCSegment > BLCoordinates::segmentRCVector [static, private] |
int BLCoordinates::currentRC [static, private] |
G4double BLCoordinates::tolerance = 1.0e-4 [static, private] |
G4double BLCoordinates::mostRecentReferenceZ = 0 [static, private] |
Referenced by getMostRecentReferenceZ(), and setMostRecentReferenceZ().
int BLCoordinates::referenceCLsegment [private] |
Referenced by BLCoordinates().
Referenced by BLCoordinates(), getCoords(), getRotation(), and updateReferenceCoordinates().