Belle II Software development
ObserverCheckFilters.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11//#include <framework/logging/Logger.h>
12#include <framework/datastore/StoreArray.h>
13
14// contains also the helper struct SpacePointInfo
15#include <tracking/dataobjects/ObserverInfo.h>
16
17#include <tracking/spacePointCreation/SpacePoint.h>
18
19// #include <memory> // unique_ptr
20
21
22namespace Belle2 {
32 public:
33
37
40
43
44
45
47 static bool initialize(const StoreArray<ObserverInfo>& aStoreArray)
48 {
50 return true;
51 }
52
54 template <class Var, class Range, typename ... types >
55 static bool initialize(Var, Range, const types& ...)
56 {
57 return true;
58 }
59
60
62 static SpacePointInfo convertSpacePoint(const SpacePoint& aSpacePoint)
63 {
64 SpacePointInfo spInfo;
65 spInfo.setPosition(aSpacePoint.getPosition());
66 spInfo.setPositionError(aSpacePoint.getPositionError());
67 spInfo.setNormalizedLocalU(aSpacePoint.getNormalizedLocalU());
68 spInfo.setNormalizedLocalV(aSpacePoint.getNormalizedLocalV());
69 spInfo.setClustersAssignedU(aSpacePoint.getIfClustersAssigned().first);
70 spInfo.setClustersAssignedV(aSpacePoint.getIfClustersAssigned().second);
71 spInfo.setVxdID(aSpacePoint.getVxdID());
72 spInfo.setSensorType((int)aSpacePoint.getType());
73 spInfo.setQualityIndicator(aSpacePoint.getQualityEstimation());
74 spInfo.setIsAssigned(aSpacePoint.getAssignmentState());
75 return spInfo;
76 }
77
78
81 static void prepare(const SpacePoint& outerHit,
82 const SpacePoint& innerHit)
83 {
84 // discard all previous data
86 std::vector<SpacePointInfo> hitsinfo = { ObserverCheckFilters::convertSpacePoint(outerHit),
88 };
89 s_observerInfo.setHits(hitsinfo);
90
91
92 // collect purity for each particle attached to the hits
93 // TODO: this is code used in the ObserverCheckMCPurity check if the following code
94 /*
95 std::vector<const Belle2::SpacePoint*> hits = {&outerHit, &innerHit};
96 std::vector<Belle2::MCVXDPurityInfo> particlesFound;
97 // TODO: check if it can cope with real data (no MC particles)!
98 particlesFound = createPurityInfosVec(hits);
99 // the dominating-particle is the uppermost one:
100 auto purityPack = particlesFound.at(0).getPurity(); //that part does not look safe to be used on data!
101 s_mainMCParticleID = purityPack.first;
102 s_mainPurity = purityPack.second;
103 */
104 }
105
108 template < typename ... types >
109 static void terminate(const types& ...)
110 {
111 }
112
113
115 template < typename ... types >
116 static void collect(const types& ...)
117 {
118 // append a copy to the storearray
119 s_storeArray.appendNew(s_observerInfo);
120 }
121
122
124 template<class Var, class RangeType>
125 static void notify(const Var&,
126 typename Var::variableType fResult,
127 const RangeType& range,
128 const typename Var::argumentType&,
129 const typename Var::argumentType&,
130 const typename Var::argumentType&)
131 {
132 // create input-container for purity-check:
133 // std::vector<const Belle2::SpacePoint*> spacePoints = {&outerHit, &centerHit, &innerHit};
134
135 generalNotify<Var, RangeType>(fResult, range);
136 }
137
138
140 template<class Var, class RangeType>
141 static void notify(const Var&,
142 typename Var::variableType fResult,
143 const RangeType& range,
144 const typename Var::argumentType&,
145 const typename Var::argumentType&)
146 {
147 // // create input-container for purity-check:
148 // std::vector<const Belle2::SpacePoint*> spacePoints = {&outerHit, &innerHit};
149
150 generalNotify<Var, RangeType>(fResult, range);
151 }
152
153
154
155
156
157
158 protected:
159
163 template<class Var, class RangeType>
164 static void generalNotify(typename Var::variableType fResult,
165 const RangeType& range)
166 {
167 // store the data retrieved:
168 FilterInfo info(Var().name(), double(fResult), range.contains(fResult), true);
170 }
171 };
172
173
175}
176
177
helper class to store the information for a Filter
Definition: FilterInfo.h:20
this observer searches logs the response for each of SelectionVariables used in the filters If the po...
ObserverCheckFilters()
empty constructor:
static void terminate(const types &...)
static method used by the observed object to terminate the observer.
static SpacePointInfo convertSpacePoint(const SpacePoint &aSpacePoint)
convert a SpacePiont into a version that can be stored in the datastore
static void prepare(const SpacePoint &outerHit, const SpacePoint &innerHit)
static method used by the observed object to reset the stored values of the observer.
static void notify(const Var &, typename Var::variableType fResult, const RangeType &range, const typename Var::argumentType &, const typename Var::argumentType &)
notifier which finds the mcParticles attached to given pair of spacePoints and determines the puritie...
static void generalNotify(typename Var::variableType fResult, const RangeType &range)
unified part of the notifier function.
static bool initialize(Var, Range, const types &...)
this function is needed by the Filters.h but has no task for this observer
static StoreArray< ObserverInfo > s_storeArray
hold a storearray to have access to the datastore
static void notify(const Var &, typename Var::variableType fResult, const RangeType &range, const typename Var::argumentType &, const typename Var::argumentType &, const typename Var::argumentType &)
notifier which finds the mcParticles attached to given triplet of spacePoints and determines the puri...
static ObserverInfo s_observerInfo
container that stores the results calculated for a selectionVariableName, has to be static due to the...
static void collect(const types &...)
fill the storearray
static bool initialize(const StoreArray< ObserverInfo > &aStoreArray)
get a copy of a storearray
Helper class that stores the information an Observer stores: i.e.
Definition: ObserverInfo.h:24
void setHits(const std::vector< SpacePointInfo > &newHits)
sets the hits the filter has been evaluated with
Definition: ObserverInfo.h:103
void clear()
resets all member variables
Definition: ObserverInfo.h:32
void addFilterInfo(FilterInfo info)
add a new filter info:
Definition: ObserverInfo.h:86
Represents a range of arithmetic types.
Definition: Range.h:29
helper class to store the SpacePoint information as coding convention prohibits to use the SpacePoint...
void setNormalizedLocalU(double val)
setter for the normalized u coordinate
void setSensorType(int type)
setter for sensor type:
void setPosition(ROOT::Math::XYZVector v)
setter for the position.
void setQualityIndicator(double qi)
setter for the quality indicator
void setVxdID(Belle2::VxdID::baseType anId)
setter for the VxdID:
void setPositionError(ROOT::Math::XYZVector v)
setter for the uncertainty on the position
void setNormalizedLocalV(double val)
setter for the normalized v coordinate
void setClustersAssignedV(bool b)
setter for is v cluster assigned
void setIsAssigned(bool ia)
setter for is assigned
void setClustersAssignedU(bool b)
setter for is u cluster assigned
SpacePoint typically is build from 1 PXDCluster or 1-2 SVDClusters.
Definition: SpacePoint.h:42
float getQualityEstimation() const
Getter for the quality of this SpacePoint.
Definition: SpacePoint.h:191
double getNormalizedLocalV() const
Return normalized local coordinates of the cluster in v (0 <= posV <= 1).
Definition: SpacePoint.h:154
std::pair< bool, bool > getIfClustersAssigned() const
Returns, if u(v)-coordinate is based on cluster information.
Definition: SpacePoint.h:162
VxdID getVxdID() const
Return the VxdID of the sensor on which the the cluster of the SpacePoint lives.
Definition: SpacePoint.h:148
const B2Vector3D & getPositionError() const
return the hitErrors in sigma of the global position
Definition: SpacePoint.h:141
bool getAssignmentState() const
Getter for status of assignment to a track.
Definition: SpacePoint.h:185
Belle2::VXD::SensorInfoBase::SensorType getType() const
Return SensorType (PXD, SVD, ...) on which the SpacePoint lives.
Definition: SpacePoint.h:145
const B2Vector3D & getPosition() const
return the position vector in global coordinates
Definition: SpacePoint.h:138
double getNormalizedLocalU() const
Return normalized local coordinates of the cluster in u (0 <= posU <= 1).
Definition: SpacePoint.h:151
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
B2Vector3D outerHit(0, 0, 0)
testing out of range behavior
Abstract base class for different kinds of events.