CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CombinedKinematicConstraintT.h
Go to the documentation of this file.
1 #ifndef CombinedKinematicConstraintT_H
2 #define CombinedKinematicConstraintT_H
3 
8 
9 #if defined( __GXX_EXPERIMENTAL_CXX0X__)
10 // this is generic: to be moved elsewhere
11 #include<tuple>
12 #include<functional>
13 #include<algorithm>
14 #include<cassert>
15 
16 // run time iteration
17 template<class TupleType, size_t N>
18 struct do_iterate
19 {
20  template<typename F>
21  static void call(TupleType& t, F f)
22  {
23  f(std::get<N-1>(t));
24  do_iterate<TupleType, N-1>::call(t,f);
25  }
26  template<typename F>
27  static void call(TupleType const & t, F f)
28  {
29  f(std::get<N-1>(t));
30  do_iterate<TupleType, N-1>::call(t,f);
31  }
32 
33 
34 };
35 
36 template<class TupleType>
37 struct do_iterate<TupleType, 0>
38 {
39  template<typename F>
40  static void call(TupleType&, F)
41  {}
42  template<typename F>
43  static void call(TupleType const &, F)
44  {}
45 };
46 
47 template<class TupleType, typename F>
48 void iterate_tuple(TupleType& t, F f)
49 {
51 }
52 
53 template<class TupleType, typename F>
54 void iterate_tuple(TupleType const& t, F f)
55 {
57 }
58 
59 
60 namespace combinedConstraintHelpers {
61 
62  // a bit less generic
64  struct totDim {
65  typedef typename std::tuple_element<N-1,TupleType>::type Elem;
66  enum { nDim = Elem::nDim + totDim<TupleType,N-1>::nDim};
67  };
68 
69  template<class TupleType>
70  struct totDim<TupleType, 0> {
71  enum { nDim=0};
72  };
73 
74  template<typename T>
75  void sum2(T& x, T y) { x+=y;}
76 
77  // mind: iteration is backward...
78  template<int DIM>
79  struct Place {
80  int offset;
81  Place() : offset(DIM) {}
82  ~Place() {
83  assert(offset==DIM || offset==0);
84  }
85  };
86 
87  template<int DIM>
88  struct PlaceValue : public Place<DIM> {
89  PlaceValue(ROOT::Math::SVector<double, DIM> & iret) : ret(iret){}
90  ROOT::Math::SVector<double, DIM> & ret;
91  template<typename C>
92  void operator()(C const & cs) {
93  this->offset -= C::nDim;
94  ret.Place_at(cs.value(),this->offset);
95  }
96  };
97 
98  template<int DIM, int NTRK>
99  struct PlaceParDer : public Place<DIM> {
100  PlaceParDer(ROOT::Math::SMatrix<double, DIM, 7*NTRK> & iret) : ret(iret){}
101  ROOT::Math::SMatrix<double, DIM, 7*NTRK> & ret;
102  template<typename C>
103  void operator()(C const & cs) {
104  this->offset -= C::nDim;
105  ret.Place_at(cs.parametersDerivative(),this->offset,0);
106  }
107  };
108 
109  template<int DIM>
110  struct PlacePosDer : public Place<DIM> {
111  PlacePosDer(ROOT::Math::SMatrix<double, DIM, 3> & iret) : ret(iret){}
112  ROOT::Math::SMatrix<double, DIM, 3> & ret;
113  template<typename C>
114  void operator()(C const & cs) {
115  this->offset -= C::nDim;
116  ret.Place_at(cs.positionDerivative(),this->offset,0);
117  }
118  };
119 
120 
121 }
122 
123 
135 // maybe a variadic template will be better
136 template< class TupleType, int NTRK >
137 class CombinedKinematicConstraintT : public MultiTrackKinematicConstraintT<NTRK, combinedConstraintHelpers::totDim<TupleType>::nDim>{
138 
139  // need compile time assert on NTRK
140 public:
143  typedef typename super::valueType valueType;
144  typedef typename super::parametersDerivativeType parametersDerivativeType;
145  typedef typename super::positionDerivativeType positionDerivativeType;
146 
147  typedef TupleType Constraints;
148 
149  //FIXME
150  enum {DIM = super::nDim};
151 
152 
153 public:
154  CombinedKinematicConstraintT(Constraints const & iconstraints): constraints(constraints){
155  }
156 
157  // initialize the constraint so it can precompute common qualtities to the three next call
158  virtual void init(const std::vector<KinematicState>& states,
159  const GlobalPoint& point, const GlobalVector& mf) {
160  iterate_tuple(constraints,
161  std::bind(&base::init,std::placeholders::_1,std::ref(states),std::ref(point), std::ref(mf)));
162  }
163 
164 
165 private:
171  void fillValue() const{
172  combinedConstraintHelpers::PlaceValue<DIM> helper(super::vl());
173  iterate_tuple(constraints,std::ref(helper));
174  }
175 
181  void fillParametersDerivative() const{
182  combinedConstraintHelpers::PlaceParDer<DIM,NTRK> helper(super::jac_d());
183  iterate_tuple(constraints,std::ref(helper));
184  }
185 
191  void fillPositionDerivative() const{
192  combinedConstraintHelpers::PlacePosDer<DIM> helper(super::jac_e());
193  iterate_tuple(constraints,std::ref(helper));
194  }
195 
196 public:
200  virtual int numberOfEquations() const {
201  int tot=0;
202  iterate_tuple(constraints,std::bind(combinedConstraintHelpers::sum2<int>,std::ref(tot),
203  std::bind(&base::numberOfEquations,std::placeholders::_1)
204  )
205  );
206  return tot;
207  }
208 
209  virtual CombinedKinematicConstraintT * clone() const
210  {
211  return new CombinedKinematicConstraintT(*this);
212  }
213 
214 private:
215  Constraints constraints;
216 
217 };
218 
219 #endif // __GXX_EXPERIMENTAL_CXX0X__
220 
221 #endif
tuple base
Main Program
Definition: newFWLiteAna.py:92
type
Definition: HCALResponse.h:21
ROOT::Math::SMatrix< double, DIM, 3 > positionDerivativeType
auto_ptr< ClusterSequence > cs
int init
Definition: HydjetWrapper.h:63
ROOT::Math::SMatrix< double, DIM, 7 *NTRK > parametersDerivativeType
virtual int numberOfEquations() const =0
double f[11][100]
unsigned int offset(bool)
ROOT::Math::SVector< double, DIM > valueType
virtual void fillPositionDerivative() const =0
#define N
Definition: blowfish.cc:9
#define DIM
virtual MultiTrackKinematicConstraintBaseT * clone() const =0
virtual void fillParametersDerivative() const =0
virtual void init(const std::vector< KinematicState > &states, const GlobalPoint &point, const GlobalVector &mf)=0
Definition: DDAxes.h:10
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
virtual void fillValue() const =0
long double T
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5