CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BasicMultiTrajectoryState.cc
Go to the documentation of this file.
2 
5 
6 using namespace SurfaceSideDefinition;
7 
8 BasicMultiTrajectoryState::BasicMultiTrajectoryState( const std::vector<TSOS>& tsvec) :
9  theCombinedStateUp2Date( false)
10 {
11  theStates.reserve(tsvec.size());
12  for (std::vector<TSOS>::const_iterator i=tsvec.begin(); i!=tsvec.end(); i++) {
13  if (!i->isValid()) {
14  throw cms::Exception("LogicError") << "MultiTrajectoryState constructed with invalid state";
15  }
16  if (i->hasError() != tsvec.front().hasError()) {
17  throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with and without errors";
18  }
19  if ( &i->surface() != &tsvec.front().surface()) {
20  throw cms::Exception("LogicError") << "MultiTrajectoryState mixes states with different surfaces";
21  }
22  if ( i->surfaceSide() != tsvec.front().surfaceSide()) {
23  throw cms::Exception("LogicError")
24  << "MultiTrajectoryState mixes states defined before and after material";
25  }
26  if ( i->localParameters().pzSign()*tsvec.front().localParameters().pzSign()<0. ) {
27  throw cms::Exception("LogicError")
28  << "MultiTrajectoryState mixes states with different signs of local p_z";
29  }
30  if ( i==tsvec.begin() ) {
31  // only accept planes!!
32  const BoundPlane* bp = dynamic_cast<const BoundPlane*>(&i->surface());
33  if ( bp==0 )
34  throw cms::Exception("LogicError") << "MultiTrajectoryState constructed on cylinder";
35  }
36  theStates.push_back( *i);
37  }
38 }
39 
41 {
42  if (theCombinedStateUp2Date) return;
43 
46 
47 }
48 
50 
51  if (theStates.empty()) {
52  edm::LogError("BasicMultiTrajectoryState")
53  << "Asking for weight of empty MultiTrajectoryState, returning zero!";
54  return 0.;
55  }
56 
57  double sumw = 0.;
58  for (std::vector<TSOS>::const_iterator it = theStates.begin(); it != theStates.end(); it++) {
59  sumw += it->weight();
60  }
61  return sumw;
62 }
63 
64 
66 
67  if (theStates.empty()) {
68  edm::LogError("BasicMultiTrajectoryState") << "Trying to rescale errors of empty MultiTrajectoryState!";
69  return;
70  }
71 
72  for (std::vector<TSOS>::iterator it = theStates.begin(); it != theStates.end(); it++) {
73  it->rescaleError(factor);
74  }
76 }
77 
78 const MagneticField*
80 {
81  //
82  // Magnetic field should be identical in all components:
83  // avoid forcing the combination of states and take value from 1st component!
84  //
85  if (theStates.empty()) {
86  edm::LogError("BasicMultiTrajectoryState")
87  << "Asking for magneticField of empty MultiTrajectoryState, returning null pointer!";
88  return 0;
89  }
90  return theStates.front().magneticField();
91 }
92 
95 {
96  //
97  // SurfaceSide should be identical in all components:
98  // avoid forcing the combination of states and take value from 1st component!
99  //
100  if (theStates.empty()) {
101  edm::LogError("BasicMultiTrajectoryState")
102  << "Asking for magneticField of empty MultiTrajectoryState, returning atCenterOfSurface!";
103  return atCenterOfSurface;
104  }
105  return theStates.front().surfaceSide();
106 }
107 
108 void
111  const Surface& aSurface,
112  const MagneticField* field,
113  const SurfaceSide side)
114 {
115  throw cms::Exception("LogicError",
116  "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, Surface, ...) called even if canUpdateLocalParameters() is false");
117 }
118 
119 void
122  const LocalTrajectoryError& err,
123  const Surface& aSurface,
124  const MagneticField* field,
125  const SurfaceSide side,
126  double weight)
127 {
128  throw cms::Exception("LogicError",
129  "BasicMultiTrajectoryState::update(LocalTrajectoryParameters, LocalTrajectoryError, ...) called even if canUpdateLocalParameters() is false");
130 }
131 
int i
Definition: DBlmapReader.cc:9
TrajectoryStateOnSurface combine(const std::vector< TrajectoryStateOnSurface > &tsos) const
const MagneticField * magneticField() const
virtual void update(const LocalTrajectoryParameters &p, const Surface &aSurface, const MagneticField *field, const SurfaceSide side)
MultiTrajectoryStateCombiner theCombiner
virtual SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.