Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
C
cbmroot
Manage
Activity
Members
Labels
Plan
Wiki
Redmine
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Container Registry
Model registry
Operate
Environments
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Computing
cbmroot
Commits
d6fa2131
Commit
d6fa2131
authored
3 years ago
by
Dominik Smith
Browse files
Options
Downloads
Patches
Plain Diff
Added time difference to closest MC event per seed to CbmSeedFinderQa.
parent
187b28be
No related branches found
Branches containing commit
No related tags found
Tags containing commit
1 merge request
!741
SeedFinderQaUpgrade
Changes
2
Hide whitespace changes
Inline
Side-by-side
Showing
2 changed files
reco/eventbuilder/digis/CbmSeedFinderQa.cxx
+35
-4
35 additions, 4 deletions
reco/eventbuilder/digis/CbmSeedFinderQa.cxx
reco/eventbuilder/digis/CbmSeedFinderQa.h
+3
-0
3 additions, 0 deletions
reco/eventbuilder/digis/CbmSeedFinderQa.h
with
38 additions
and
4 deletions
reco/eventbuilder/digis/CbmSeedFinderQa.cxx
+
35
−
4
View file @
d6fa2131
...
...
@@ -39,9 +39,13 @@ CbmSeedFinderQa::CbmSeedFinderQa() : fOutFolder("SeedFinderQA", "Seed finder QA"
fhCorrectVsFoundNoNoise
=
new
TH2I
(
"fhCorrectVsFoundNoNoise"
,
"Correct digis [pct] vs. Found digis [pct], no noise; Correct; Found "
,
110
,
-
5.
,
105.
,
110
,
-
5.
,
105.
);
fhTimeOffset
=
new
TH1F
(
"fhTimeOffset"
,
"tSeed - tMCEvent [ns]"
,
20
,
-
5
,
5
);
fhTimeOffset
=
new
TH1F
(
"fhTimeOffsetMatched"
,
"tSeed - tMCMatched [ns]"
,
20
,
-
5
,
5
);
fhTimeOffset
->
SetCanExtend
(
TH1
::
kAllAxes
);
fhTimeOffsetClosest
=
new
TH1F
(
"fhTimeOffsetClosest"
,
"tSeed - tMCClosest [ns]"
,
20
,
-
5
,
5
);
fhTimeOffsetClosest
->
SetCanExtend
(
TH1
::
kAllAxes
);
histFolder
->
Add
(
fhCorrectDigiRatio
);
histFolder
->
Add
(
fhCorrectDigiRatioNoNoise
);
histFolder
->
Add
(
fhNoiseDigiRatio
);
...
...
@@ -49,6 +53,7 @@ CbmSeedFinderQa::CbmSeedFinderQa() : fOutFolder("SeedFinderQA", "Seed finder QA"
histFolder
->
Add
(
fhCorrectVsFound
);
histFolder
->
Add
(
fhCorrectVsFoundNoNoise
);
histFolder
->
Add
(
fhTimeOffset
);
histFolder
->
Add
(
fhTimeOffsetClosest
);
histFolder
->
Add
(
fhLinkedMCEventsPerTrigger
);
histFolder
->
Add
(
fhLinkedTriggersPerMCEvent
);
histFolder
->
Add
(
fhMatchedTriggersPerMCEvent
);
...
...
@@ -67,6 +72,7 @@ CbmSeedFinderQa::~CbmSeedFinderQa()
delete
fhCorrectVsFound
;
delete
fhCorrectVsFoundNoNoise
;
delete
fhTimeOffset
;
delete
fhTimeOffsetClosest
;
delete
fhLinkedMCEventsPerTrigger
;
delete
fhLinkedTriggersPerMCEvent
;
delete
fhMatchedTriggersPerMCEvent
;
...
...
@@ -82,12 +88,17 @@ void CbmSeedFinderQa::Init()
fEventList
=
(
CbmMCEventList
*
)
FairRootManager
::
Instance
()
->
GetObject
(
"MCEventList."
);
}
void
CbmSeedFinderQa
::
ResetPerTsStorage
()
{
fvEventMatchesPerTs
.
clear
();
}
void
CbmSeedFinderQa
::
ResetPerTsStorage
()
{
fvEventMatchesPerTs
.
clear
();
fvSeedTimesPerTs
.
clear
();
}
void
CbmSeedFinderQa
::
FillQaSeedInfo
(
const
int32_t
WinStart
,
const
int32_t
WinEnd
,
const
std
::
vector
<
CbmMatch
>*
vDigiMatch
,
const
double
seedTime
)
{
fvSeedTimesFull
.
push_back
(
seedTime
);
fvSeedTimesPerTs
.
push_back
(
seedTime
);
int32_t
digiCount
=
0
;
int32_t
noiseDigiCount
=
0
;
...
...
@@ -178,6 +189,8 @@ void CbmSeedFinderQa::FillQaSeedInfo(const int32_t WinStart, const int32_t WinEn
void
CbmSeedFinderQa
::
FillQaMCInfo
()
{
const
uint32_t
nEvents
=
fEventList
->
GetNofEvents
();
if
(
nEvents
==
0
)
{
return
;
}
std
::
vector
<
uint32_t
>
vLinkedTriggersPerMCEvent
;
std
::
vector
<
uint32_t
>
vMatchedTriggersPerMCEvent
;
vLinkedTriggersPerMCEvent
.
resize
(
nEvents
,
0
);
...
...
@@ -200,6 +213,21 @@ void CbmSeedFinderQa::FillQaMCInfo()
for
(
const
auto
&
value
:
vMatchedTriggersPerMCEvent
)
{
fhMatchedTriggersPerMCEvent
->
Fill
(
value
);
}
// get sorted vector of MC event times
std
::
vector
<
double
>
vMCEventTimes
;
for
(
uint32_t
iEvent
=
0
;
iEvent
<
nEvents
;
iEvent
++
)
{
vMCEventTimes
.
push_back
(
fEventList
->
GetEventTimeByIndex
(
iEvent
));
}
std
::
sort
(
std
::
begin
(
vMCEventTimes
),
std
::
end
(
vMCEventTimes
));
//find closest MC event for each seed (assumes both arrays are sorted in time)
auto
minElem
=
vMCEventTimes
.
begin
();
for
(
const
auto
&
seedTime
:
fvSeedTimesPerTs
)
{
auto
comp
=
[
&
,
seedTime
](
double
val1
,
double
val2
)
{
return
fabs
(
seedTime
-
val1
)
<
fabs
(
seedTime
-
val2
);
};
minElem
=
std
::
min_element
(
minElem
,
vMCEventTimes
.
end
(),
comp
);
fhTimeOffsetClosest
->
Fill
(
seedTime
-
*
minElem
);
}
}
void
CbmSeedFinderQa
::
FillHistos
()
...
...
@@ -278,12 +306,15 @@ void CbmSeedFinderQa::WriteHistos()
fhTimeOffset
->
DrawCopy
(
"colz"
,
""
);
fCanv
->
cd
(
8
);
fh
LinkedMCEventsPerTrigger
->
DrawCopy
(
"colz"
,
""
);
fh
TimeOffsetClosest
->
DrawCopy
(
"colz"
,
""
);
fCanv
->
cd
(
9
);
fhLinked
TriggersPerMCEvent
->
DrawCopy
(
"colz"
,
""
);
fhLinked
MCEventsPerTrigger
->
DrawCopy
(
"colz"
,
""
);
fCanv
->
cd
(
10
);
fhLinkedTriggersPerMCEvent
->
DrawCopy
(
"colz"
,
""
);
fCanv
->
cd
(
11
);
fhMatchedTriggersPerMCEvent
->
DrawCopy
(
"colz"
,
""
);
FairSink
*
sink
=
FairRootManager
::
Instance
()
->
GetSink
();
...
...
This diff is collapsed.
Click to expand it.
reco/eventbuilder/digis/CbmSeedFinderQa.h
+
3
−
0
View file @
d6fa2131
...
...
@@ -63,6 +63,8 @@ public:
/** @brief Matches that link constructed event seeds to MC events, current timeslice only. */
std
::
vector
<
CbmMatch
>
fvEventMatchesPerTs
;
/** @brief Vector of event seeds, current TS only. */
std
::
vector
<
double
>
fvSeedTimesPerTs
;
/** @brief Gather QA Information.
* @params WinStart Starting position of seed window.
...
...
@@ -103,6 +105,7 @@ private:
TH2I
*
fhCorrectVsFound
=
nullptr
;
/// correct digis per event vs found digis per event
TH2I
*
fhCorrectVsFoundNoNoise
=
nullptr
;
/// correct digis per event vs found digis per event, disregarding noise
TH1F
*
fhTimeOffset
=
nullptr
;
/// difference between true event time and seed time
TH1F
*
fhTimeOffsetClosest
=
nullptr
;
/// difference between seed time and closest MC event time
CbmQaCanvas
*
fCanv
;
///summary canvas
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
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!
Save comment
Cancel
Please
register
or
sign in
to comment