Public Member Functions | Private Member Functions | Private Attributes | List of all members
geo::PlaneGeo Class Reference

Geometry information for a single readout plane. More...

#include "/cvmfs/nova-development.opensciencegrid.org/novasoft/releases/N20-09-17/GeometryObjects/PlaneGeo.h"

Public Member Functions

 PlaneGeo (std::vector< const TGeoNode * > &n, unsigned int depth)
 Construct a representation of a single plane of the detector. More...
 
 ~PlaneGeo ()
 
unsigned int Ncells () const
 Number of cells in this plane. More...
 
const CellGeoCell (int icell) const
 
View_t View () const
 Which coordinate does this plane measure. More...
 
void LocalToWorld (const double *plane, double *world) const
 Transform point from local plane frame to world frame. More...
 
void LocalToWorldVect (const double *plane, double *world) const
 Transform direction vector from local to world. More...
 
void WorldToLocal (const double *world, double *plane) const
 Transform point from world frame to local plane frame. More...
 
void WorldToLocalVect (const double *world, double *plane) const
 Transform direction vector from world to local. More...
 
void TranslatePlane (double dx, double dy=0, double dz=0, double phi=0, double theta=0, double psi=0)
 
void TranslatePlane (TGeoCombiTrans *comb)
 
void RestoreOriginal ()
 Restore geometry to original. More...
 

Private Member Functions

void FindCells (std::vector< const TGeoNode * > &n, unsigned int depth)
 
void MakeCell (std::vector< const TGeoNode * > &n, unsigned int depth)
 

Private Attributes

TGeoHMatrix * fGeoMatrix
 Plane to world transform. More...
 
TGeoHMatrix * fGeoMatrixOriginal
 Original plane to world transform. More...
 
bool fIsOriginal
 Is fGeoMatrix the original transform? More...
 
Readout_t fReadout
 Which end is the readout on? More...
 
View_t fView
 Does this plane measure X or Y? More...
 
std::vector< CellGeo * > fCell
 List of cells in this plane. More...
 

Detailed Description

Geometry information for a single readout plane.

Definition at line 36 of file PlaneGeo.h.

Constructor & Destructor Documentation

geo::PlaneGeo::PlaneGeo ( std::vector< const TGeoNode * > &  n,
unsigned int  depth 
)

Construct a representation of a single plane of the detector.

Definition at line 42 of file PlaneGeo.cxx.

References a1, stan::math::fabs(), fCell, fGeoMatrix, fGeoMatrixOriginal, FindCells(), fIsOriginal, fReadout, fView, MECModelEnuComparisons::i, geo::kEast, geo::kTop, geo::kWest, geo::kX, geo::kY, geo::sort_hori(), geo::sort_vert(), and fillBadChanDBTables::step.

43  {
44  // build a matrix to take us from the local to the world coordinates
45  // in one step
46  fGeoMatrix = new TGeoHMatrix(*n[0]->GetMatrix());
47  for (unsigned int i=1; i<=depth; ++i) {
48  fGeoMatrix->Multiply(n[i]->GetMatrix());
49  }
50 
51  // Copy fGeoMatrix to Original
52  fGeoMatrixOriginal = new TGeoHMatrix(*fGeoMatrix);
53 
54  // fGeoMatrix is still unmodified
55  fIsOriginal = true;
56 
57  // Determine read-out orientation for the plane. By convention, the
58  // planes, cells, and fibers are constructed such that the z-axis
59  // goes along the long direction with z increasing as one moves
60  // toward the read out end. So, a good way to find the read out end
61  // is to take a small step in the +z direction which moves you
62  // toward the read out side of the detector.
63  static const double step = 10.0; // size of step in cm
64  static const double tol = 0.1*step; // testing tolerance
65  static const double a0[3] = {0,0,0};
66  static const double a1[3] = {0,0,step};
67  double b0[3], b1[3];
68  fGeoMatrix->LocalToMaster(a0,b0);
69  fGeoMatrix->LocalToMaster(a1,b1);
70  b1[0] -= b0[0]; b1[1] -= b0[1]; b1[2] -= b0[2];
71  if (fabs(b1[0]-0.0)<tol && fabs(b1[1]-step)<tol) {
72  fView = kX;
73  fReadout = kTop;
74  }
75  else if (fabs(b1[0]-step)<tol && fabs(b1[1]-0.)<tol) {
76  fView = kY;
77  fReadout = kWest;
78  }
79  else if (fabs(b1[0]+step)<tol && fabs(b1[1]-0.)<tol) {
80  fView = kY;
81  fReadout = kEast;
82  }
83  else {
84  throw cet::exception("BAD_PLANE_VIEW")
85  << "Could not determine view for plane\n"
86  << __FILE__ << ":" << __LINE__ << "\n";
87  }
88  this->FindCells(n, depth);
89 
90  // Make sure horizontals are sorted top to bottom and verticals are
91  // sorted east to west
92  if (fView == kX) {
93  std::sort(fCell.begin(), fCell.end(), sort_vert);
94  }
95  else if (fView == kY) {
96  std::sort(fCell.begin(), fCell.end(), sort_hori);
97  }
98  else abort();
99  }
fvar< T > fabs(const fvar< T > &x)
Definition: fabs.hpp:15
void FindCells(std::vector< const TGeoNode * > &n, unsigned int depth)
Definition: PlaneGeo.cxx:115
Vertical planes which measure X.
Definition: PlaneGeo.h:28
Vertical modules readout on the top.
Definition: PlaneGeo.h:21
::xsd::cxx::tree::exception< char > exception
Definition: Database.h:225
View_t fView
Does this plane measure X or Y?
Definition: PlaneGeo.h:88
Horizontal modules readout on the west.
Definition: PlaneGeo.h:22
Horizontal planes which measure Y.
Definition: PlaneGeo.h:29
bool fIsOriginal
Is fGeoMatrix the original transform?
Definition: PlaneGeo.h:86
TGeoHMatrix * fGeoMatrixOriginal
Original plane to world transform.
Definition: PlaneGeo.h:85
static bool sort_vert(const CellGeo *c1, const CellGeo *c2)
Definition: PlaneGeo.cxx:32
Readout_t fReadout
Which end is the readout on?
Definition: PlaneGeo.h:87
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
TH1F * a1
Definition: f2_nu.C:476
static bool sort_hori(const CellGeo *c1, const CellGeo *c2)
Definition: PlaneGeo.cxx:25
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89
Horizontal modules readout on the east.
Definition: PlaneGeo.h:23
geo::PlaneGeo::~PlaneGeo ( )

Definition at line 102 of file PlaneGeo.cxx.

References plot_validation_datamc::c, fCell, fGeoMatrix, and fGeoMatrixOriginal.

103  {
104  if(fGeoMatrix) delete fGeoMatrix;
106 
107  for(size_t c = 0; c < fCell.size(); ++c){
108  if(fCell[c]) { delete fCell[c]; fCell[c] = 0; }
109  }
110  fCell.clear();
111  }
TGeoHMatrix * fGeoMatrixOriginal
Original plane to world transform.
Definition: PlaneGeo.h:85
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89

Member Function Documentation

const CellGeo* geo::PlaneGeo::Cell ( int  icell) const
inline

Return the icell'th cell in the plane. Should start from bottom and count up in case of horizontals. Should count from west to east in case of verticals

Definition at line 48 of file PlaneGeo.h.

References fCell.

Referenced by dif::DiFShowerFinder::adjustPlane(), dif::DiFShowerFinder::adjustPlane_end(), geo::GeometryTest::analyze(), zcl::SMMTriggerAna::analyze(), align::Alignment::analyze(), mcchk::CosmicAna::analyze(), mcchk::RockAna::analyze(), sn::SNSlicerAna::analyze(), calib::AssessCalib::analyze(), calib::HitEfficiency::analyze(), zcl::FmmTrackerValidation::analyze(), trk::KalmanTrackAna::analyze(), remid::ReMIdDedxRock::analyze(), remid::ReMIdDedxFD::analyze(), remid::ReMIdDedxStudies::analyze(), LightLevels::analyze(), htk::HoughTrack::analyze(), trk::CosmicTrackAna::analyze(), calib::CalibAna::ApproximateHitPos(), trk::CosmicTrackAna::beginJob(), calib::BestPathEstimates(), geo::GeometryBase::BuildMaps(), slid::ParticleIDAlg::CalcTrkHitPath(), trk::KalmanTrackMerge::CanJoinTracks(), geo::GeometryBase::CellInfo(), dt::Chunk::CellPoints(), geo::GeometryBase::CellTpos(), geo::GeometryBase::ClosestApproach(), dt::Cand::ClosestToEndCell(), geo::GeometryBase::CountCellsOnLineFast(), trk::KalmanGeoHelper::CountMissedCellsOnLine(), LightLevels::CreateHistos(), murem::TrackCleanUpAlg::DeDxInPlane(), evd::GeometryDrawer::DrawCells2D(), evd::GeometryDrawer::DrawHighlightCell(), skim::ParametersNumu::EarliestLatestHitPos(), mcchk::ShowerAnaCheck::EMMoliereRadius(), mcchk::ShowerAnaCheck::EMShowerWidth(), calib::BetheBlochFit::endJob(), dif::DiFShowerFinder::eparm(), calib::MakeAttenuationProfiles::event(), geo::LiveGeometry::FillBadBoxes(), calib::StopperThreshold::FillHist(), trk::KalmanGeoHelper::fillprivategeo(), trk::CosmicTrackAna::FillTrackHistograms(), calib::StopperThreshold::FillTree(), calib::CosmicCalib::filter(), bsf::BremShowerFilter::filter(), calib::FindBelowThresholdCalibCandidates(), trk::WindowTrackingAlg::FindEndPoint(), bsf::BremShowerFilter::findShowerByReco(), dif::DiFShowerFinder::findShowerByReco(), bpfit::DimuonFitter::FindVertexZ(), trk::CosmicTrackUtilities::FindZBoundaries(), calib::FindZBoundaries(), bpfit::DimuonFitter::FitView(), trk::WindowTrackingAlg::FitWindow(), align::SplitTracks::GeoEstimatePos(), calib::Calibrator::GetAttenCurve(), evd::GeometryDrawer::GetBox(), jmshower::RecoJMShower::GetCellDistToPoint(), jmshower::RecoJMShower::GetCellDistToTrk(), align::SplitTracks::GetCellEndpoints(), slid::ParticleIDAlg::GetCellNodePos(), jmshower::RecoJMShower::GetCellNodePos(), calib::CosmicCalib::GetChannelHists(), calib::MakeAttenuationProfiles::GetChannelHists(), nuesand::FillNueSandbox::GetdEdx(), geo::LiveGeometry::GetDetectorEdges(), skim::ParametersNumu::GetDetectorEdges(), airshower::AirSlicer::GetHitPos(), jmshower::JMShower::GetTransHits(), jmshower::RecoJMShower::GetTrkHitPath(), jmshower::RecoJMShower::GetTrkHitPos(), jmshower::RecoJMShower::GetTrkPlaneCell(), jmshower::RecoJMShower::GetTrkPlaneDistToEdge(), jmshower::RecoJMShower::GetTrkPlanePos(), calib::GetXYZD(), dt::Chunk::HitsOnLine(), nerd::NERDProng::HitToHitDistance(), me::MEFinder::HitToHitDistance(), geo::GeometryBase::IdToCell(), bsf::BremShowerFilter::inFiducial(), dif::DiFShowerFinder::inFiducial(), skim::EvaluatorNue::KeepSlice(), calib::AttenCache::LoadFromDatabase(), trk::CosmicTrackAlg::MakeTrack(), trk::WindowTrackingAlg::MakeTrack(), trk::WindowTrackingAlg::MakeViewTrack(), calib::MakeZBoundaryMap(), fuzz::ViewMatchAlg::Matching(), nerd::ViewMatchAlg::Matching(), lem::FindLEMMatches::MatchToVertex(), lem::LEM::MatchToVertex(), trk::KalmanGeoHelper::MatchTrajectoryToPlane(), zcl::SMMCluster::MissC(), dif::DiFShowerFinder::muonstub(), numue::NumuEAlg::NumuEAlg(), skim::ParametersNumu::ParametersNumu(), trk::CosmicTrackUtilities::PathLengthInCell(), calib::PathLengthInCell(), geo::plane_sort(), rb::Track::PlaneDirMap(), slid::ParticleIDAlg::PlaneHitCell(), slid::ParticleIDAlg::PlaneHitXYZ(), slid::Recluster::produce(), jmshower::JMTrackMerge::produce(), zcl::SMMCluster::produce(), photrans::PhotonTransport::produce(), cvn::RegCVNMapper::produce(), cheat::MCCheater::produce(), remid::RecoMuon::produce(), caf::CAFMaker::produce(), evd::RawDataDrawer::RawDigit2D(), jmshower::RecoJMShower::RecoShowers(), trk::RLFit::RLFit(), hough::MultiHoughT::Scrub(), earms::ElasticArmsHS::Scrub(), vdt::VertexDT::Scrub(), fuzz::FuzzyKVertex::Scrub(), trk::KalmanTrackMerge::ShiftInterpolationPoints(), trk::WindowTrackingAlg::ShortTrackExtraPlane(), trk::WindowTrackingAlg::ShortViewTrack(), murem::TrackCleanUpAlg::SortByDistFromTrack(), photrans::ImprovedTransport::StepAlongHit(), geo::GeometryTest::testCellId(), geo::GeometryTest::testCellIdFromPos(), geo::GeometryTest::testCellPos(), geo::GeometryTest::testFindCell(), calib::StopperThreshold::testPath(), geo::GeometryTest::testStepping(), geo::GeometryTest::testUniqueId(), trk::CosmicTrackUtilities::TruePathLengthInCell(), slid::NuEEnergyAlg::VertexEnergy(), lem::VertexToPlaneAndCell(), trk::KalmanTrackMerge::ViewMergeTracks(), and calib::ZBounds().

48  {
49  if (icell < 0 || icell >= (int)fCell.size()) return 0;
50  return fCell[icell]; }
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89
void geo::PlaneGeo::FindCells ( std::vector< const TGeoNode * > &  n,
unsigned int  depth 
)
private

Definition at line 115 of file PlaneGeo.cxx.

References GetName(), MECModelEnuComparisons::i, MakeCell(), nd, and registry_explorer::v.

Referenced by PlaneGeo(), and View().

117  {
118  if (strncmp(n[depth]->GetName(),"vCell",5)==0) {
119  this->MakeCell(n,depth);
120  return;
121  }
122 
123  // Explore the next layer down
124  unsigned int deeper = depth+1;
125  if (deeper>=n.size()) {
126  throw cet::exception("BAD_NODE")
127  << "Exceeded maximum depth\n"
128  << __FILE__ << ":" << __LINE__ << "\n";
129  }
130  const TGeoVolume* v = n[depth]->GetVolume();
131  int nd = v->GetNdaughters();
132  for (int i=0; i<nd; ++i) {
133  n[deeper] = v->GetNode(i);
134  this->FindCells(n, deeper);
135  }
136  }
void FindCells(std::vector< const TGeoNode * > &n, unsigned int depth)
Definition: PlaneGeo.cxx:115
::xsd::cxx::tree::exception< char > exception
Definition: Database.h:225
std::string GetName(int i)
void MakeCell(std::vector< const TGeoNode * > &n, unsigned int depth)
Definition: PlaneGeo.cxx:140
void geo::PlaneGeo::LocalToWorld ( const double *  plane,
double *  world 
) const

Transform point from local plane frame to world frame.

Definition at line 148 of file PlaneGeo.cxx.

References fGeoMatrix.

Referenced by View().

149  {
150  fGeoMatrix->LocalToMaster(plane, world);
151  }
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
void geo::PlaneGeo::LocalToWorldVect ( const double *  plane,
double *  world 
) const

Transform direction vector from local to world.

Definition at line 155 of file PlaneGeo.cxx.

References fGeoMatrix.

Referenced by View().

156  {
157  fGeoMatrix->LocalToMasterVect(plane, world);
158  }
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
void geo::PlaneGeo::MakeCell ( std::vector< const TGeoNode * > &  n,
unsigned int  depth 
)
private

Definition at line 140 of file PlaneGeo.cxx.

References fCell.

Referenced by FindCells(), and View().

142  {
143  fCell.push_back(new CellGeo(n, depth));
144  }
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89
unsigned int geo::PlaneGeo::Ncells ( ) const
inline

Number of cells in this plane.

Definition at line 43 of file PlaneGeo.h.

References fCell.

Referenced by geo::GeometryTest::analyze(), mcchk::CosmicAna::analyze(), numue::NumuEAna::analyze(), align::Alignment::analyze(), htk::HoughTrack::analyze(), chaninfo::BadChanValidate::BadChansFromDB(), chaninfo::CosmicEff::beginRun(), dprf::ChannelPlots::beginRun(), dprf::TrackPlots::beginRun(), calib::DriftResponseCalc::beginRun(), rsim::ReadoutSim::beginRun(), calib::CosmicTrends::beginRun(), chaninfo::BadChanValidate::beginSubRun(), noe::build_cell_lookup_table(), geo::GeometryBase::BuildMaps(), slid::ParticleIDAlg::CalcCellPlaneTransverseDedx(), slid::ParticleIDAlg::CalcPlaneTransverseDedx(), geo::GeometryBase::CellInfo(), geo::GeometryBase::CellTpos(), remid::ReMIdTrain::ContainedEvent(), remid::ReMIdDedx::ContainedEvent(), remid::ReMIdDedxRock::ContainedEvent(), remid::ReMIdDedxStudies::ContainedEvent(), geo::GeometryBase::CountCellsOnLineFast(), trk::KalmanGeoHelper::CountMissedCellsOnLine(), LightLevels::CreateHistos(), rsim::ReadoutSim::CreateRawDigits(), evd::GeometryDrawer::DrawCells2D(), chaninfo::CosmicEff::endJob(), calib::AttenFit::endJob(), comi::ChanOcc::endSubRun(), rsim::RecordNoiseSpectrumFile::endSubRun(), LightLevels::endSubRun(), rsim::IFPGAAlgorithm::FetchThresholds(), caf::FillHadClustVars(), trk::KalmanGeoHelper::fillprivategeo(), caf::FillSliceInfo(), caf::FillSliceVars(), evd::GeometryDrawer::GetBox(), jmshower::RecoJMShower::GetCellTransDedx(), evd::GeometryDrawer::GetDCMBoxes(), jmshower::RecoJMShower::GetTrkCPlaneCell(), jmshower::RecoJMShower::GetTrkPlaneCell(), calib::GetXYZD(), calib::HasXYAdjacents(), dt::DiscreteTracker::HighActivity(), photrans::ImprovedTransport::LoadHistos(), chaninfo::BadChanList::LoadNewData(), lem::FindLEMMatches::MatchToVertex(), lem::LEM::MatchToVertex(), trk::KalmanGeoHelper::MatchTrajectoryToPlane(), skim::ParametersNumu::ParametersNumu(), slid::ParticleIDAlg::PlaneHitCell(), chaninfo::BadChanList::preBeginEvent(), trk::KalmanTrack::produce(), numusand::FillSandbox::produce(), lem::AlignLibToVtx::produce(), cvn::RegCVNMapper::produce(), rsim::ReadoutSim::produce(), evd::RawDataDrawer::RawDigit2D(), lem::GenFromLib::readNext(), chaninfo::BadChanList::SetRandomBadChannels(), geo::GeometryTest::testCellId(), geo::GeometryTest::testCellIdFromPos(), geo::GeometryTest::testCellPos(), lem::VertexToPlaneAndCell(), dt::View::View(), and calib::AttenuationFit::writeResults().

43 { return fCell.size(); }
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89
void geo::PlaneGeo::RestoreOriginal ( )

Restore geometry to original.

Restore original matrix

Definition at line 216 of file PlaneGeo.cxx.

References fGeoMatrix, fGeoMatrixOriginal, and fIsOriginal.

Referenced by View().

217  {
218 
219  /// Restore original matrix
220  fGeoMatrix = new TGeoHMatrix(*fGeoMatrixOriginal);
221 
222  // fGeoMatrix is back to original
223  fIsOriginal = true;
224 
225  }
bool fIsOriginal
Is fGeoMatrix the original transform?
Definition: PlaneGeo.h:86
TGeoHMatrix * fGeoMatrixOriginal
Original plane to world transform.
Definition: PlaneGeo.h:85
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
void geo::PlaneGeo::TranslatePlane ( double  dx,
double  dy = 0,
double  dz = 0,
double  phi = 0,
double  theta = 0,
double  psi = 0 
)

Translate plane and all cells in plane by dx, dy and dz Original transform will be stored in a separate TGeoHMatrix

phi, theta and psi are euler angles

dx, dy and dz are the translations we want

Definition at line 179 of file PlaneGeo.cxx.

Referenced by align::Alignment::analyze(), and View().

181  {
182 
183  /// phi, theta and psi are euler angles
184  TGeoRotation *rot = new TGeoRotation("alignmentRot",phi,theta,psi);
185  /// dx, dy and dz are the translations we want
186  TGeoTranslation *trans = new TGeoTranslation("",dx,dy,dz);
187  TGeoCombiTrans *comb = new TGeoCombiTrans(*trans,*rot);
188 
189  this->TranslatePlane(comb);
190 
191  }
double dy[NP][NC]
double dx[NP][NC]
double dz[NP][NC]
void TranslatePlane(double dx, double dy=0, double dz=0, double phi=0, double theta=0, double psi=0)
Definition: PlaneGeo.cxx:179
void geo::PlaneGeo::TranslatePlane ( TGeoCombiTrans *  comb)

Translate plane using a previously defined TGeoCombiTrans Original transform will be stored in a separate TGeoHMatrix

Restart from original matrix

Apply shift

Now apply shift to all cells

Definition at line 195 of file PlaneGeo.cxx.

References fCell, fGeoMatrix, fGeoMatrixOriginal, fIsOriginal, and MECModelEnuComparisons::i.

196  {
197 
198  /// Restart from original matrix
199  fGeoMatrix = new TGeoHMatrix(*fGeoMatrixOriginal);
200 
201  /// Apply shift
202  fGeoMatrix->MultiplyLeft(comb->MakeClone());
203 
204  /// Now apply shift to all cells
205  for(unsigned int i=0; i<fCell.size(); i++){
206  fCell.at(i)->TranslateCell(comb);
207  }
208 
209  // fGeoMatrix is now modified
210  fIsOriginal = false;
211 
212  }
bool fIsOriginal
Is fGeoMatrix the original transform?
Definition: PlaneGeo.h:86
TGeoHMatrix * fGeoMatrixOriginal
Original plane to world transform.
Definition: PlaneGeo.h:85
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
std::vector< CellGeo * > fCell
List of cells in this plane.
Definition: PlaneGeo.h:89
View_t geo::PlaneGeo::View ( ) const
inline

Which coordinate does this plane measure.

Definition at line 53 of file PlaneGeo.h.

References dx, dy, dz, FindCells(), fView, LocalToWorld(), LocalToWorldVect(), MakeCell(), NDAPDHVSetting::plane, RestoreOriginal(), chisquared::theta, TranslatePlane(), WorldToLocal(), and WorldToLocalVect().

Referenced by dprf::ChannelPlots::analyze(), mcchk::DetAna::analyze(), slicer::S4DParamCalc::analyze(), dprf::TrackPlots::analyze(), align::Alignment::analyze(), calib::HitEfficiency::analyze(), remid::ReMIdDedxRock::analyze(), remid::ReMIdDedxFD::analyze(), remid::ReMIdDedxStudies::analyze(), calib::CosmicTrends::analyze(), comi::Leana::AnaNHit(), dprf::ChannelPlots::beginRun(), calib::BestPathEstimates(), slid::ParticleIDAlg::CalcAsymIneria(), geo::GeometryBase::CellInfo(), dt::Chunk::CellPoints(), geo::GeometryBase::CellTpos(), murem::TrackCleanUpAlg::CleanUpTrack(), geo::GeometryBase::ClosestApproach(), geo::GeometryBase::CountCellsOnLineFast(), trk::KalmanGeoHelper::CountMissedCellsOnLine(), util::CountXY(), murem::TrackCleanUpAlg::DeDxInPlane(), lem::DefaultVertex(), evd::GeometryDrawer::DrawCells2D(), evd::GeometryDrawer::DrawHighlightCell(), calib::FindOutliers::endRun(), comi::ChanOcc::endSubRun(), util::EventBox(), calib::CosmicTrends::FillHistograms(), bpfit::BreakPoint::FillHitList(), calib::FindBelowThresholdCalibCandidates(), bpfit::BreakPoint::FitTracks(), bpfit::DimuonFitter::FitView(), calib::Calibrator::GetAttenCurve(), evd::GeometryDrawer::GetBox(), jmshower::RecoJMShower::GetCellDistToPoint(), jmshower::RecoJMShower::GetCellDistToTrk(), slid::ParticleIDAlg::GetCellNodePos(), jmshower::RecoJMShower::GetCellNodePos(), evd::GeometryDrawer::GetDCMBoxes(), airshower::AirSlicer::GetHitPos(), jmshower::RecoJMShower::GetTrkHitPos(), jmshower::RecoJMShower::GetTrkPlaneCell(), calib::ThresholdAna::GetView(), calib::CosmicCalib::GetView(), calib::MakeAttenuationProfiles::GetView(), calib::AttenFit::GetView(), calib::AttenuationFit::GetView(), calib::AttenCache::GetView(), calib::GetXYZD(), rb::HitList::HitList(), calib::AttenCache::LoadFromDatabase(), evtsum::EventSummary::MakeOutput(), lem::FindLEMMatches::MatchToVertex(), lem::LEM::MatchToVertex(), trk::KalmanGeoHelper::MatchTrajectoryToPlane(), trk::KalmanGeoHelper::MatchTrajectoryToPlaneInView(), trk::CosmicTrackUtilities::PathLengthInCell(), calib::PathLengthInCell(), slid::ParticleIDAlg::PlaneHitCell(), trk::KalmanTrack::produce(), lem::MakeLibrary::produce(), cvn::RegCVNMapper::produce(), cheat::MCCheater::produce(), remid::RecoMuon::produce(), evd::RawDataDrawer::RawDigit2D(), jmshower::RecoJMShower::RecoShowers(), calib::AttenCache::RepresentativePlane(), trk::RLFit::RLFit(), trk::WindowTrackingAlg::SetTrackEndPoints(), photrans::ImprovedTransport::StepAlongHit(), geo::GeometryTest::testCellPos(), and geo::GeometryTest::testFindPlanes().

53 { return fView; }
View_t fView
Does this plane measure X or Y?
Definition: PlaneGeo.h:88
void geo::PlaneGeo::WorldToLocal ( const double *  world,
double *  plane 
) const

Transform point from world frame to local plane frame.

Definition at line 162 of file PlaneGeo.cxx.

References fGeoMatrix.

Referenced by View().

163  {
164  fGeoMatrix->MasterToLocal(world, plane);
165  }
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84
void geo::PlaneGeo::WorldToLocalVect ( const double *  world,
double *  plane 
) const

Transform direction vector from world to local.

Convert a vector from world frame to the local plane frame

Parameters
world: 3-D array. Vector in world coordinates; input.
plane: 3-D array. Vector in plane coordinates; plane.

Definition at line 172 of file PlaneGeo.cxx.

References fGeoMatrix.

Referenced by View().

173  {
174  fGeoMatrix->MasterToLocalVect(world,plane);
175  }
TGeoHMatrix * fGeoMatrix
Plane to world transform.
Definition: PlaneGeo.h:84

Member Data Documentation

std::vector<CellGeo*> geo::PlaneGeo::fCell
private

List of cells in this plane.

Definition at line 89 of file PlaneGeo.h.

Referenced by Cell(), MakeCell(), Ncells(), PlaneGeo(), TranslatePlane(), and ~PlaneGeo().

TGeoHMatrix* geo::PlaneGeo::fGeoMatrix
private

Plane to world transform.

Definition at line 84 of file PlaneGeo.h.

Referenced by LocalToWorld(), LocalToWorldVect(), PlaneGeo(), RestoreOriginal(), TranslatePlane(), WorldToLocal(), WorldToLocalVect(), and ~PlaneGeo().

TGeoHMatrix* geo::PlaneGeo::fGeoMatrixOriginal
private

Original plane to world transform.

Definition at line 85 of file PlaneGeo.h.

Referenced by PlaneGeo(), RestoreOriginal(), TranslatePlane(), and ~PlaneGeo().

bool geo::PlaneGeo::fIsOriginal
private

Is fGeoMatrix the original transform?

Definition at line 86 of file PlaneGeo.h.

Referenced by PlaneGeo(), RestoreOriginal(), and TranslatePlane().

Readout_t geo::PlaneGeo::fReadout
private

Which end is the readout on?

Definition at line 87 of file PlaneGeo.h.

Referenced by PlaneGeo().

View_t geo::PlaneGeo::fView
private

Does this plane measure X or Y?

Definition at line 88 of file PlaneGeo.h.

Referenced by PlaneGeo(), and View().


The documentation for this class was generated from the following files: