Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Sign in
Toggle navigation
Menu
Open sidebar
Corryvreckan
Corryvreckan
Commits
c54fdfbd
Commit
c54fdfbd
authored
Jun 24, 2021
by
Simon Spannagel
Browse files
Merge branch 'checkTracking' into 'master'
Adding a module `AnalysisTracks` See merge request
!444
parents
77e2265c
949e1635
Pipeline
#2777418
passed with stages
in 17 minutes and 29 seconds
Changes
4
Pipelines
8
Hide whitespace changes
Inline
Side-by-side
src/modules/AnalysisTracks/AnalysisTracks.cpp
0 → 100644
View file @
c54fdfbd
/**
* @file
* @brief Implementation of module AnalysisTracks
*
* @copyright Copyright (c) 2020 CERN and the Corryvreckan authors.
* This software is distributed under the terms of the MIT License, copied verbatim in the file "LICENSE.md".
* In applying this license, CERN does not waive the privileges and immunities granted to it by virtue of its status as an
* Intergovernmental Organization or submit itself to any jurisdiction.
*/
#include
"AnalysisTracks.h"
using
namespace
corryvreckan
;
AnalysisTracks
::
AnalysisTracks
(
Configuration
&
config
,
std
::
vector
<
std
::
shared_ptr
<
Detector
>>
detectors
)
:
Module
(
config
,
std
::
move
(
detectors
))
{}
void
AnalysisTracks
::
initialize
()
{
std
::
string
title
;
for
(
auto
&
detector
:
get_detectors
())
{
if
(
detector
->
isAuxiliary
())
{
continue
;
}
LOG
(
DEBUG
)
<<
"Initialise for detector "
+
detector
->
getName
();
TDirectory
*
directory
=
getROOTDirectory
();
TDirectory
*
local_directory
=
directory
->
mkdir
(
detector
->
getName
().
c_str
());
if
(
local_directory
==
nullptr
)
{
throw
RuntimeError
(
"Cannot create or access local ROOT directory for module "
+
this
->
getUniqueName
());
}
local_directory
->
cd
();
title
=
"Distance between tracks ; distance [mm]; entries"
;
_distance_between_tracks_
[
detector
->
getName
().
c_str
()]
=
new
TH1F
(
"distance_tracks"
,
title
.
c_str
(),
1000
,
0
,
10
);
title
=
"Tracks with same hits ; # tracks with same hit; entries"
;
_tracks_per_hit_
[
detector
->
getName
().
c_str
()]
=
new
TH1F
(
"number_tracks_with same hit"
,
title
.
c_str
(),
15
,
0
,
15
);
title
=
"Cluster vs tracks ; # tracks; #clusters"
;
clusters_vs_tracks_
[
detector
->
getName
().
c_str
()]
=
new
TH2F
(
"clusters_vs_tracks"
,
title
.
c_str
(),
25
,
0
,
25
,
200
,
0
,
200
);
}
}
StatusCode
AnalysisTracks
::
run
(
const
std
::
shared_ptr
<
Clipboard
>&
clipboard
)
{
auto
tracks
=
clipboard
->
getData
<
Track
>
();
if
(
!
tracks
.
size
())
{
return
StatusCode
::
Success
;
}
for
(
auto
d
:
get_detectors
())
{
if
(
d
->
isAuxiliary
())
{
continue
;
}
clusters_vs_tracks_
.
at
(
d
->
getName
())
->
Fill
(
static_cast
<
double
>
(
tracks
.
size
()),
static_cast
<
double
>
(
clipboard
->
getData
<
Cluster
>
(
d
->
getName
()).
size
()));
}
// Loop over all tracks and get clusters assigned to tracks as well as the intersections
std
::
map
<
std
::
string
,
std
::
map
<
std
::
pair
<
double
,
double
>
,
int
>>
track_clusters
;
std
::
map
<
std
::
string
,
std
::
vector
<
XYZPoint
>>
intersects
;
// local coordinates
for
(
auto
&
track
:
tracks
)
{
for
(
auto
d
:
get_detectors
())
{
if
(
d
->
isAuxiliary
())
{
continue
;
}
intersects
[
d
->
getName
()].
push_back
(
d
->
globalToLocal
(
track
->
getState
(
d
->
getName
())));
if
(
d
->
isDUT
()
||
track
->
getClusterFromDetector
(
d
->
getName
())
==
nullptr
)
{
continue
;
}
auto
c
=
track
->
getClusterFromDetector
(
d
->
getName
());
track_clusters
[
d
->
getName
()][
std
::
make_pair
<
double
,
double
>
(
c
->
column
(),
c
->
row
())]
+=
1
;
}
}
// Now fill the histos
for
(
auto
const
&
intersect
:
intersects
)
{
auto
key
=
intersect
.
first
;
auto
val
=
intersect
.
second
;
for
(
uint
i
=
0
;
i
<
val
.
size
();
++
i
)
{
auto
j
=
i
+
1
;
while
(
j
<
val
.
size
())
{
auto
p
=
val
.
at
(
i
)
-
val
.
at
(
j
);
_distance_between_tracks_
.
at
(
key
)
->
Fill
((
p
.
Mag2
()));
j
++
;
}
}
}
for
(
auto
const
&
track_cluster
:
track_clusters
)
{
auto
key
=
track_cluster
.
first
;
for
(
auto
const
&
v
:
track_cluster
.
second
)
{
_tracks_per_hit_
.
at
(
key
)
->
Fill
(
v
.
second
);
}
}
// Return value telling analysis to keep running
return
StatusCode
::
Success
;
}
src/modules/AnalysisTracks/AnalysisTracks.h
0 → 100644
View file @
c54fdfbd
/**
* @file
* @brief Definition of module AnalysisTracks
*
* @copyright Copyright (c) 2020 CERN and the Corryvreckan authors.
* This software is distributed under the terms of the MIT License, copied verbatim in the file "LICENSE.md".
* In applying this license, CERN does not waive the privileges and immunities granted to it by virtue of its status as an
* Intergovernmental Organization or submit itself to any jurisdiction.
*/
#ifndef ANALYSISTRACKS_H
#define ANALYSISTRACKS_H 1
#include
<TCanvas.h>
#include
<TH1F.h>
#include
<TH2F.h>
#include
<iostream>
#include
"core/module/Module.hpp"
#include
"objects/Cluster.hpp"
#include
"objects/Pixel.hpp"
#include
"objects/Track.hpp"
namespace
corryvreckan
{
/** @ingroup Modules
* @brief Module to do function
*
* More detailed explanation of module
*/
class
AnalysisTracks
:
public
Module
{
public:
/**
* @brief Constructor for this unique module
* @param config Configuration object for this module as retrieved from the steering file
* @param detectors Vector of pointers to the detectors
*/
AnalysisTracks
(
Configuration
&
config
,
std
::
vector
<
std
::
shared_ptr
<
Detector
>>
detectors
);
~
AnalysisTracks
()
{}
/**
* @brief [Initialise this module]
*/
void
initialize
()
override
;
/**
* @brief [Run the function of this module]
*/
StatusCode
run
(
const
std
::
shared_ptr
<
Clipboard
>&
clipboard
)
override
;
private:
std
::
map
<
std
::
string
,
TH1F
*>
_distance_between_tracks_
{};
std
::
map
<
std
::
string
,
TH1F
*>
_tracks_per_hit_
{};
std
::
map
<
std
::
string
,
TH2F
*>
clusters_vs_tracks_
{};
};
}
// namespace corryvreckan
#endif // ANALYSISTRACKS_H
src/modules/AnalysisTracks/CMakeLists.txt
0 → 100644
View file @
c54fdfbd
CORRYVRECKAN_ENABLE_DEFAULT
(
ON
)
# Define module and return the generated name as MODULE_NAME
CORRYVRECKAN_GLOBAL_MODULE
(
MODULE_NAME
)
# Add source files to library
CORRYVRECKAN_MODULE_SOURCES
(
${
MODULE_NAME
}
AnalysisTracks.cpp
# ADD SOURCE FILES HERE...
)
# Provide standard install target
CORRYVRECKAN_MODULE_INSTALL
(
${
MODULE_NAME
}
)
src/modules/AnalysisTracks/README.md
0 → 100644
View file @
c54fdfbd
# AnalysisTracks
**Maintainer**
: Lennart Huth (lennart.huth@desy.de)
**Module Type**
:
*GLOBAL*
**Status**
: Functional
### Description
Study reconstructed tracks in each event by calculating the distance between
tracks in an event. In addition cluster assignments to multiple tracks is
plotted. Finally the number of tracks is correlated with the number of tracks
reconstructed.
### Parameters
No parameters are used from the configuration file.
### Plots produced
For each detector the following plots are produced:
*
2D histogram of number of tracks versus clusters
*
1D histogram of absolute distance between tracks
*
1D histogram of number of times a cluster is used in a track
### Usage
```
toml
[AnalysisTracks]
```
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment