CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
SurveyInputCSCfromPins Class Reference

#include <SurveyInputCSCfromPins.h>

Inheritance diagram for SurveyInputCSCfromPins:
SurveyInputBase edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 Read ideal tracker geometry from DB. More...
 
 SurveyInputCSCfromPins (const edm::ParameterSet &)
 
- Public Member Functions inherited from SurveyInputBase
void beginJob () override
 Read data from input. More...
 
 ~SurveyInputBase () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void errors (double a, double b, bool missing1, bool missing2, double &dx_dx, double &dy_dy, double &dz_dz, double &phix_phix, double &phiz_phiz, double &dy_phix)
 
void fillAllRecords (Alignable *ali)
 
void orient (align::LocalVector LC1, align::LocalVector LC2, double a, double b, double &T, double &dx, double &dy, double &dz, double &PhX, double &PhZ)
 

Private Attributes

double m_errorX
 
double m_errorY
 
double m_errorZ
 
double m_missingErrorAngle
 
double m_missingErrorTranslation
 
std::string m_pinPositions
 
std::string m_rootFile
 
double m_stationErrorPhiX
 
double m_stationErrorPhiY
 
double m_stationErrorPhiZ
 
double m_stationErrorX
 
double m_stationErrorY
 
double m_stationErrorZ
 
bool m_verbose
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from SurveyInputBase
static void addComponent (Alignable *)
 
static Alignabledetector ()
 Get alignable detector as read from input. More...
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
- Protected Attributes inherited from SurveyInputBase
bool theFirstEvent
 

Detailed Description

Class to read ideal tracker from DB.

$Date: Fri Jun 29 09:20:52 CEST 2007

Revision
1.2
Author
Dmitry Yakorev

Definition at line 16 of file SurveyInputCSCfromPins.h.

Constructor & Destructor Documentation

SurveyInputCSCfromPins::SurveyInputCSCfromPins ( const edm::ParameterSet cfg)

Definition at line 25 of file SurveyInputCSCfromPins.cc.

26  : m_pinPositions(cfg.getParameter<std::string>("pinPositions"))
27  , m_rootFile(cfg.getParameter<std::string>("rootFile"))
28  , m_verbose(cfg.getParameter<bool>("verbose"))
29  , m_errorX(cfg.getParameter<double>("errorX"))
30  , m_errorY(cfg.getParameter<double>("errorY"))
31  , m_errorZ(cfg.getParameter<double>("errorZ"))
32  , m_missingErrorTranslation(cfg.getParameter<double>("missingErrorTranslation"))
33  , m_missingErrorAngle(cfg.getParameter<double>("missingErrorAngle"))
34  , m_stationErrorX(cfg.getParameter<double>("stationErrorX"))
35  , m_stationErrorY(cfg.getParameter<double>("stationErrorY"))
36  , m_stationErrorZ(cfg.getParameter<double>("stationErrorZ"))
37  , m_stationErrorPhiX(cfg.getParameter<double>("stationErrorPhiX"))
38  , m_stationErrorPhiY(cfg.getParameter<double>("stationErrorPhiY"))
39  , m_stationErrorPhiZ(cfg.getParameter<double>("stationErrorPhiZ"))
40 {}
T getParameter(std::string const &) const

Member Function Documentation

void SurveyInputCSCfromPins::analyze ( const edm::Event ,
const edm::EventSetup iSetup 
)
overridevirtual

Read ideal tracker geometry from DB.

Implements SurveyInputBase.

Definition at line 112 of file SurveyInputCSCfromPins.cc.

References a, SurveyInputBase::addComponent(), AlignableNavigator::alignableFromDetId(), align::AlignableMuon, b, AlignableMuon::CSCEndcaps(), PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, relativeConstraints::error, errors(), edmOneToOneComparison::file1, fillAllRecords(), edm::EventSetup::get(), h, mps_fire::i, cuy::ii, recoMuon::in, m_missingErrorAngle, m_pinPositions, m_rootFile, m_verbose, Alignable::mother(), Alignable::move(), orient(), Alignable::rotateAroundLocalX(), Alignable::rotateAroundLocalZ(), Alignable::setSurvey(), Alignable::surface(), Alignable::survey(), SurveyInputBase::theFirstEvent, AlignableSurface::toGlobal(), AlignableSurface::toLocal(), PV3DBase< T, PVType, FrameType >::x(), globals_cff::x1, globals_cff::x2, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

113 {
114  if (theFirstEvent) {
115 
116  edm::LogInfo("SurveyInputCSCfromPins") << "***************ENTERING INITIALIZATION******************" << " \n";
117 
118  std::ifstream in;
119  in.open(m_pinPositions.c_str());
120 
121  Double_t x1, y1, z1, x2, y2, z2, a, b, tot=0.0, maxErr=0.0, h, s1, dx, dy, dz, PhX, PhZ, T;
122 
123  int ID1, ID2, ID3, ID4, ID5, i=1, ii=0;
124 
125  TFile *file1 = new TFile(m_rootFile.c_str(),"recreate");
126  TTree *tree1 = new TTree("tree1","alignment pins");
127 
128  if (m_verbose) {
129 
130  tree1->Branch("displacement_x_pin1_cm", &x1, "x1/D");
131  tree1->Branch("displacement_y_pin1_cm", &y1, "y1/D");
132  tree1->Branch("displacement_z_pin1_cm", &z1, "z1/D");
133  tree1->Branch("displacement_x_pin2_cm", &x2, "x2/D");
134  tree1->Branch("displacement_y_pin2_cm", &y2, "y2/D");
135  tree1->Branch("displacement_z_pin2_cm", &z2, "z2/D");
136  tree1->Branch("error_vector_length_cm",&h, "h/D");
137  tree1->Branch("stretch_diff_cm",&s1, "s1/D");
138  tree1->Branch("stretch_factor",&T, "T/D");
139  tree1->Branch("chamber_displacement_x_cm",&dx, "dx/D");
140  tree1->Branch("chamber_displacement_y_cm",&dy, "dy/D");
141  tree1->Branch("chamber_displacement_z_cm",&dz, "dz/D");
142  tree1->Branch("chamber_rotation_x_rad",&PhX, "PhX/D");
143  tree1->Branch("chamber_rotation_z_rad",&PhZ, "PhZ/D");
144  }
145 
146  edm::ESHandle<DTGeometry> dtGeometry;
147  edm::ESHandle<CSCGeometry> cscGeometry;
148  iSetup.get<MuonGeometryRecord>().get( dtGeometry );
149  iSetup.get<MuonGeometryRecord>().get( cscGeometry );
150 
151  AlignableMuon* theAlignableMuon = new AlignableMuon( &(*dtGeometry) , &(*cscGeometry) );
152  AlignableNavigator* theAlignableNavigator = new AlignableNavigator( theAlignableMuon );
153 
154  const auto& theEndcaps = theAlignableMuon->CSCEndcaps();
155 
156  for (const auto& aliiter: theEndcaps) {
157  addComponent(aliiter);
158  }
159 
160 
161  while (in.good())
162  {
163  bool missing1 = false;
164  bool missing2 = false;
165 
166  in >> ID1 >> ID2 >> ID3 >> ID4 >> ID5 >> x1 >> y1 >> z1 >> x2 >> y2 >> z2 >> a >> b;
167 
168  if (fabs(x1 - 1000.) < 1e5 && fabs(y1 - 1000.) < 1e5 && fabs(z1 - 1000.) < 1e5) {
169  missing1 = true;
170  x1 = x2;
171  y1 = y2;
172  z1 = z2;
173  }
174 
175  if (fabs(x2 - 1000.) < 1e5 && fabs(y2 - 1000.) < 1e5 && fabs(z2 - 1000.) < 1e5) {
176  missing2 = true;
177  x2 = x1;
178  y2 = y1;
179  z2 = z1;
180  }
181 
182  x1=x1/10.0;
183  y1=y1/10.0;
184  z1=z1/10.0;
185  x2=x2/10.0;
186  y2=y2/10.0;
187  z2=z2/10.0;
188 
189  CSCDetId layerID(ID1, ID2, ID3, ID4, 1);
190 
191  // We cannot use chamber ID (when ID5=0), because AlignableNavigator gives the error (aliDet and aliDetUnit are undefined for chambers)
192 
193 
194  Alignable* theAlignable1 = theAlignableNavigator->alignableFromDetId( layerID );
195  Alignable* chamberAli=theAlignable1->mother();
196 
197  LocalVector LC1 = chamberAli->surface().toLocal(GlobalVector(x1,y1,z1));
198  LocalVector LC2 = chamberAli->surface().toLocal(GlobalVector(x2,y2,z2));
199 
200  orient(LC1, LC2, a, b, T, dx, dy, dz, PhX, PhZ);
201 
202  GlobalPoint PG1 = chamberAli->surface().toGlobal(LocalPoint(LC1.x(), LC1.y() + a, LC1.z() + b));
203  chamberAli->surface().toGlobal(LocalPoint(LC2.x(), LC2.y() - a, LC2.z() + b));
204 
205 
206  LocalVector lvector( dx, dy, dz);
207  GlobalVector gvector = ( chamberAli->surface()).toGlobal( lvector );
208  chamberAli->move( gvector );
209 
210  chamberAli->rotateAroundLocalX( PhX );
211  chamberAli->rotateAroundLocalZ( PhZ );
212 
213  double dx_dx, dy_dy, dz_dz, phix_phix, phiz_phiz, dy_phix;
214  errors(a, b, missing1, missing2, dx_dx, dy_dy, dz_dz, phix_phix, phiz_phiz, dy_phix);
215  align::ErrorMatrix error = ROOT::Math::SMatrixIdentity();
216  error(0,0) = dx_dx;
217  error(1,1) = dy_dy;
218  error(2,2) = dz_dz;
219  error(3,3) = phix_phix;
220  error(4,4) = m_missingErrorAngle;
221  error(5,5) = phiz_phiz;
222  error(1,3) = dy_phix;
223  error(3,1) = dy_phix; // just in case
224 
225  chamberAli->setSurvey( new SurveyDet(chamberAli->surface(), error) );
226 
227  if (m_verbose) {
228 
229  edm::LogInfo("SurveyInputCSCfromPins") << " survey information = " << chamberAli->survey() << " \n";
230 
231  LocalPoint LP1n = chamberAli->surface().toLocal(PG1);
232 
233  LocalPoint hiP(LP1n.x(), LP1n.y() - a*T, LP1n.z() - b);
234 
235  h=hiP.mag();
236  s1=LP1n.y() - a;
237 
238  if (h>maxErr) { maxErr=h;
239 
240  ii=i;
241  }
242 
243  edm::LogInfo("SurveyInputCSCfromPins") << " \n"
244  << "i " << i++ << " " << ID1 << " " << ID2 << " " << ID3 << " " << ID4 << " " << ID5 << " error " << h << " \n"
245  <<" x1 " << x1 << " y1 " << y1 << " z1 " << z1 << " x2 " << x2 << " y2 " << y2 << " z2 " << z2 << " \n"
246  << " error " << h << " S1 " << s1 << " \n"
247  <<" dx " << dx << " dy " << dy << " dz " << dz << " PhX " << PhX << " PhZ " << PhZ << " \n";
248 
249  tot+=h;
250 
251  tree1->Fill();
252  }
253 
254  }
255 
256  in.close();
257 
258  if (m_verbose) {
259 
260  file1->Write();
261  edm::LogInfo("SurveyInputCSCfromPins") << " Total error " << tot << " Max Error " << maxErr << " N " << ii << " \n";
262  }
263 
264  file1->Close();
265 
266  for (const auto& aliiter: theEndcaps) {
267  fillAllRecords(aliiter);
268  }
269 
270  delete theAlignableMuon;
271  delete theAlignableNavigator;
272 
273  edm::LogInfo("SurveyInputCSCfromPins") << "*************END INITIALIZATION***************" << " \n";
274 
275  theFirstEvent = false;
276  }
277 }
void errors(double a, double b, bool missing1, bool missing2, double &dx_dx, double &dy_dy, double &dz_dz, double &phix_phix, double &phiz_phiz, double &dy_phix)
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
virtual void rotateAroundLocalZ(Scalar radians)
Rotation around local z-axis.
Definition: Alignable.cc:234
const SurveyDet * survey() const
Return survey info.
Definition: Alignable.h:225
T y() const
Definition: PV3DBase.h:63
virtual void move(const GlobalVector &displacement)=0
Movement with respect to the global reference frame.
static void addComponent(Alignable *)
align::Alignables CSCEndcaps()
align::RotationType toLocal(const align::RotationType &) const
Return in local frame a rotation given in global frame.
void fillAllRecords(Alignable *ali)
T z() const
Definition: PV3DBase.h:64
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:135
ii
Definition: cuy.py:590
void setSurvey(const SurveyDet *)
Set survey info.
Definition: Alignable.cc:343
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
T get() const
Definition: EventSetup.h:62
virtual void rotateAroundLocalX(Scalar radians)
Rotation around local x-axis.
Definition: Alignable.cc:182
align::GlobalPoints toGlobal(const align::LocalPoints &) const
Return in global coord given a set of local points.
math::Error< 6 >::type ErrorMatrix
Definition: Definitions.h:39
long double T
Constructor of the full muon geometry.
Definition: AlignableMuon.h:37
T x() const
Definition: PV3DBase.h:62
void orient(align::LocalVector LC1, align::LocalVector LC2, double a, double b, double &T, double &dx, double &dy, double &dz, double &PhX, double &PhZ)
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:94
Global3DVector GlobalVector
Definition: GlobalVector.h:10
void SurveyInputCSCfromPins::errors ( double  a,
double  b,
bool  missing1,
bool  missing2,
double &  dx_dx,
double &  dy_dy,
double &  dz_dz,
double &  phix_phix,
double &  phiz_phiz,
double &  dy_phix 
)
private

Definition at line 67 of file SurveyInputCSCfromPins.cc.

References PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, mps_fire::i, m_errorX, m_errorY, m_errorZ, m_missingErrorTranslation, and orient().

Referenced by analyze(), and fillAllRecords().

67  {
68  dx_dx = 0.;
69  dy_dy = 0.;
70  dz_dz = 0.;
71  phix_phix = 0.;
72  phiz_phiz = 0.;
73  dy_phix = 0.;
74 
75  const double trials = 10000.; // two significant digits
76 
77  for (int i = 0; i < trials; i++) {
78  LocalVector LC1, LC2;
79  if (missing1) {
80  LC1 = LocalVector(gRandom->Gaus(0., m_missingErrorTranslation), gRandom->Gaus(0., m_missingErrorTranslation), gRandom->Gaus(0., m_missingErrorTranslation));
81  }
82  else {
83  LC1 = LocalVector(gRandom->Gaus(0., m_errorX), gRandom->Gaus(0., m_errorY), gRandom->Gaus(0., m_errorZ));
84  }
85 
86  if (missing2) {
87  LC2 = LocalVector(gRandom->Gaus(0., m_missingErrorTranslation), gRandom->Gaus(0., m_missingErrorTranslation), gRandom->Gaus(0., m_missingErrorTranslation));
88  }
89  else {
90  LC2 = LocalVector(gRandom->Gaus(0., m_errorX), gRandom->Gaus(0., m_errorY), gRandom->Gaus(0., m_errorZ));
91  }
92 
93  double dx, dy, dz, PhX, PhZ, T;
94  orient(LC1, LC2, a, b, T, dx, dy, dz, PhX, PhZ);
95 
96  dx_dx += dx * dx;
97  dy_dy += dy * dy;
98  dz_dz += dz * dz;
99  phix_phix += PhX * PhX;
100  phiz_phiz += PhZ * PhZ;
101  dy_phix += dy * PhX; // the only non-zero off-diagonal element
102  }
103 
104  dx_dx /= trials;
105  dy_dy /= trials;
106  dz_dz /= trials;
107  phix_phix /= trials;
108  phiz_phiz /= trials;
109  dy_phix /= trials;
110 }
Local3DVector LocalVector
Definition: LocalVector.h:12
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
long double T
void orient(align::LocalVector LC1, align::LocalVector LC2, double a, double b, double &T, double &dx, double &dy, double &dz, double &PhX, double &PhZ)
void SurveyInputCSCfromPins::fillAllRecords ( Alignable ali)
private

Definition at line 280 of file SurveyInputCSCfromPins.cc.

References a, funct::abs(), b, Alignable::components(), DEFINE_FWK_MODULE, MillePedeFileConverter_cfg::e, relativeConstraints::error, errors(), Alignable::geomDetId(), m_missingErrorAngle, m_missingErrorTranslation, m_stationErrorPhiX, m_stationErrorPhiY, m_stationErrorPhiZ, m_stationErrorX, m_stationErrorY, m_stationErrorZ, Alignable::setSurvey(), Alignable::surface(), and Alignable::survey().

Referenced by analyze().

280  {
281  if (ali->survey() == nullptr) {
282 
283  AlignableCSCChamber *ali_AlignableCSCChamber = dynamic_cast<AlignableCSCChamber*>(ali);
284  AlignableCSCStation *ali_AlignableCSCStation = dynamic_cast<AlignableCSCStation*>(ali);
285 
286  if (ali_AlignableCSCChamber != nullptr) {
287  CSCDetId detid(ali->geomDetId());
288  if (abs(detid.station()) == 1 && (detid.ring() == 1 || detid.ring() == 4)) {
289  align::ErrorMatrix error = ROOT::Math::SMatrixIdentity();
293  error(3,3) = m_missingErrorAngle;
294  error(4,4) = m_missingErrorAngle;
295  error(5,5) = m_missingErrorAngle;
296 
297  ali->setSurvey(new SurveyDet(ali->surface(), error));
298  }
299  else {
300  double a = 100.;
301  double b = -9.4034;
302  if (abs(detid.station()) == 1 && detid.ring() == 2) a = -90.260;
303  else if (abs(detid.station()) == 1 && detid.ring() == 3) a = -85.205;
304  else if (abs(detid.station()) == 2 && detid.ring() == 1) a = -97.855;
305  else if (abs(detid.station()) == 2 && detid.ring() == 2) a = -164.555;
306  else if (abs(detid.station()) == 3 && detid.ring() == 1) a = -87.870;
307  else if (abs(detid.station()) == 3 && detid.ring() == 2) a = -164.555;
308  else if (abs(detid.station()) == 4 && detid.ring() == 1) a = -77.890;
309 
310  double dx_dx, dy_dy, dz_dz, phix_phix, phiz_phiz, dy_phix;
311  errors(a, b, true, true, dx_dx, dy_dy, dz_dz, phix_phix, phiz_phiz, dy_phix);
312  align::ErrorMatrix error = ROOT::Math::SMatrixIdentity();
313  error(0,0) = dx_dx;
314  error(1,1) = dy_dy;
315  error(2,2) = dz_dz;
316  error(3,3) = phix_phix;
317  error(4,4) = m_missingErrorAngle;
318  error(5,5) = phiz_phiz;
319  error(1,3) = dy_phix;
320  error(3,1) = dy_phix; // just in case
321 
322  ali->setSurvey(new SurveyDet(ali->surface(), error));
323  }
324  }
325 
326  else if (ali_AlignableCSCStation != nullptr) {
327  align::ErrorMatrix error = ROOT::Math::SMatrixIdentity();
328  error(0,0) = m_stationErrorX;
329  error(1,1) = m_stationErrorY;
330  error(2,2) = m_stationErrorZ;
331  error(3,3) = m_stationErrorPhiX;
332  error(4,4) = m_stationErrorPhiY;
333  error(5,5) = m_stationErrorPhiZ;
334 
335  ali->setSurvey(new SurveyDet(ali->surface(), error));
336  }
337 
338  else {
339  align::ErrorMatrix error = ROOT::Math::SMatrixIdentity();
340  ali->setSurvey(new SurveyDet(ali->surface(), error*(1e-10)));
341  }
342  }
343 
344  for (const auto& iter: ali->components()) {
345  fillAllRecords(iter);
346  }
347 }
void errors(double a, double b, bool missing1, bool missing2, double &dx_dx, double &dy_dy, double &dz_dz, double &phix_phix, double &phiz_phiz, double &dy_phix)
const SurveyDet * survey() const
Return survey info.
Definition: Alignable.h:225
virtual const Alignables & components() const =0
Return vector of all direct components.
void fillAllRecords(Alignable *ali)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:135
void setSurvey(const SurveyDet *)
Set survey info.
Definition: Alignable.cc:343
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
math::Error< 6 >::type ErrorMatrix
Definition: Definitions.h:39
const DetId & geomDetId() const
Definition: Alignable.h:186
A muon CSC Chamber( an AlignableDet )
void SurveyInputCSCfromPins::orient ( align::LocalVector  LC1,
align::LocalVector  LC2,
double  a,
double  b,
double &  T,
double &  dx,
double &  dy,
double &  dz,
double &  PhX,
double &  PhZ 
)
private

Definition at line 42 of file SurveyInputCSCfromPins.cc.

References a, b, PV3DBase< T, PVType, FrameType >::mag(), SQR, mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by analyze(), and errors().

42  {
43  double cosPhX, sinPhX, cosPhZ, sinPhZ;
44 
45  LocalPoint LP1(LC1.x(), LC1.y() + a, LC1.z() + b);
46  LocalPoint LP2(LC2.x(), LC2.y() - a, LC2.z() + b);
47 
48  LocalPoint P((LP1.x() - LP2.x())/(2.*a), (LP1.y() - LP2.y())/(2.*a), (LP1.z() - LP2.z())/(2.*a));
49  LocalPoint Pp((LP1.x() + LP2.x())/(2.), (LP1.y() + LP2.y())/(2.), (LP1.z() + LP2.z())/(2.));
50 
51  T=P.mag();
52 
53  sinPhX=P.z()/T;
54  cosPhX=sqrt(1-SQR(sinPhX));
55  cosPhZ=P.y()/(T*cosPhX);
56  sinPhZ=-P.x()/(T*cosPhZ);
57 
58  PhX=atan2(sinPhX,cosPhX);
59 
60  PhZ=atan2(sinPhZ,cosPhZ);
61 
62  dx=Pp.x()-sinPhZ*sinPhX*b;
63  dy=Pp.y()+cosPhZ*sinPhX*b;
64  dz=Pp.z()-cosPhX*b;
65 }
T y() const
Definition: PV3DBase.h:63
T mag() const
Definition: PV3DBase.h:67
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
double b
Definition: hdecay.h:120
std::pair< OmniClusterRef, TrackingParticleRef > P
double a
Definition: hdecay.h:121
long double T
T x() const
Definition: PV3DBase.h:62
#define SQR(x)

Member Data Documentation

double SurveyInputCSCfromPins::m_errorX
private

Definition at line 36 of file SurveyInputCSCfromPins.h.

Referenced by errors().

double SurveyInputCSCfromPins::m_errorY
private

Definition at line 36 of file SurveyInputCSCfromPins.h.

Referenced by errors().

double SurveyInputCSCfromPins::m_errorZ
private

Definition at line 36 of file SurveyInputCSCfromPins.h.

Referenced by errors().

double SurveyInputCSCfromPins::m_missingErrorAngle
private

Definition at line 37 of file SurveyInputCSCfromPins.h.

Referenced by analyze(), and fillAllRecords().

double SurveyInputCSCfromPins::m_missingErrorTranslation
private

Definition at line 37 of file SurveyInputCSCfromPins.h.

Referenced by errors(), and fillAllRecords().

std::string SurveyInputCSCfromPins::m_pinPositions
private

Definition at line 33 of file SurveyInputCSCfromPins.h.

Referenced by analyze().

std::string SurveyInputCSCfromPins::m_rootFile
private

Definition at line 34 of file SurveyInputCSCfromPins.h.

Referenced by analyze().

double SurveyInputCSCfromPins::m_stationErrorPhiX
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

double SurveyInputCSCfromPins::m_stationErrorPhiY
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

double SurveyInputCSCfromPins::m_stationErrorPhiZ
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

double SurveyInputCSCfromPins::m_stationErrorX
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

double SurveyInputCSCfromPins::m_stationErrorY
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

double SurveyInputCSCfromPins::m_stationErrorZ
private

Definition at line 38 of file SurveyInputCSCfromPins.h.

Referenced by fillAllRecords().

bool SurveyInputCSCfromPins::m_verbose
private

Definition at line 35 of file SurveyInputCSCfromPins.h.

Referenced by analyze().