Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Update & fix to SeedingHitSet #36819

Merged
merged 2 commits into from
Feb 1, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
40 changes: 24 additions & 16 deletions RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h
Original file line number Diff line number Diff line change
Expand Up @@ -17,37 +17,45 @@ class SeedingHitSet {

SeedingHitSet() = default;

SeedingHitSet(ConstRecHitPointer one, ConstRecHitPointer two) : theRecHits({one, two}) { setSize(); }
SeedingHitSet(ConstRecHitPointer one, ConstRecHitPointer two) : SeedingHitSet({one, two}) {}

SeedingHitSet(ConstRecHitPointer one, ConstRecHitPointer two, ConstRecHitPointer three)
: theRecHits({one, two, three}) {
setSize();
}
: SeedingHitSet({one, two, three}) {}

SeedingHitSet(ConstRecHitPointer one, ConstRecHitPointer two, ConstRecHitPointer three, ConstRecHitPointer four)
: theRecHits({one, two, three, four}) {
setSize();
: SeedingHitSet({one, two, three, four}) {}

SeedingHitSet(const std::vector<ConstRecHitPointer> &hits) {
if (hits.size() >= 2) {
auto end = std::find(hits.begin(), hits.end(), nullPtr());
auto size = std::distance(hits.begin(), end);
if (size >= 2) {
theRecHits.reserve(size);
std::copy(hits.begin(), end, std::back_inserter(theRecHits));
}
}
}

SeedingHitSet(const std::vector<ConstRecHitPointer> &hits) : theRecHits(hits) { setSize(); }
SeedingHitSet(const std::initializer_list<ConstRecHitPointer> &hits) {
if (hits.size() >= 2) {
auto end = std::find(hits.begin(), hits.end(), nullPtr());
auto size = std::distance(hits.begin(), end);
if (size >= 2) {
theRecHits.reserve(size);
std::copy(hits.begin(), end, std::back_inserter(theRecHits));
}
}
}

ConstRecHitPointer const *data() const { return theRecHits.data(); }

unsigned int size() const { return theSize; }
unsigned int size() const { return theRecHits.size(); }

ConstRecHitPointer get(unsigned int i) const { return theRecHits[i]; }
ConstRecHitPointer operator[](unsigned int i) const { return theRecHits[i]; }

protected:
std::vector<ConstRecHitPointer> theRecHits;
unsigned int theSize = 0;

void setSize() {
theSize = 0;
while (theRecHits[++theSize] and theSize < theRecHits.size())
;
theSize = theSize > 1 ? theSize : 0;
}
};

#endif
5 changes: 5 additions & 0 deletions RecoTracker/TkSeedingLayers/test/BuildFile.xml
Original file line number Diff line number Diff line change
Expand Up @@ -5,3 +5,8 @@
<library file="*.cc" name="RecoTrackerTkSeedingLayersTests">
<flags EDM_PLUGIN="1"/>
</library>

<bin file="test_catch2_*.cpp" name="test_catch2_RECOTRACKER_TKSEEDINGLAYERS">
<use name="DataFormats/TrackerRecHit2D"/>
<use name="catch2"/>
</bin>
74 changes: 74 additions & 0 deletions RecoTracker/TkSeedingLayers/test/test_catch2_SeedingHitSet.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,74 @@
#include "RecoTracker/TkSeedingLayers/interface/SeedingHitSet.h"
#include "DataFormats/TrackerRecHit2D/interface/TrackerSingleRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
#include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
#include "DataFormats/TrackerRecHit2D/interface/MTDTrackingRecHit.h"
#include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"

#include "catch.hpp"

#include <iostream>

TEST_CASE("SeedingHitSet testing constructors", "[SeedingHitSet]") {
using Hit = SeedingHitSet::RecHit;
using HitPointer = SeedingHitSet::ConstRecHitPointer;

HitPointer hitOne = new SiPixelRecHit();
HitPointer hitTwo = new MTDTrackingRecHit();
HitPointer hitThree = new SiStripRecHit1D();
HitPointer hitFour = new SiStripRecHit2D();
HitPointer hitFive = new Phase2TrackerRecHit1D();

HitPointer HIT_NULL = nullptr;

SECTION("Check vector constructor") {
std::vector<HitPointer> containerOne = {hitOne, hitTwo, hitThree, hitFour, hitFive};
std::vector<HitPointer> containerTwo = {hitOne, hitTwo, hitThree, HIT_NULL, hitFive};
std::vector<HitPointer> containerThree = {hitOne, HIT_NULL, hitThree, hitFour, hitFive};

SeedingHitSet vecOne(containerOne);
REQUIRE(vecOne.size() == containerOne.size());
SeedingHitSet vecTwo(containerTwo);
REQUIRE(vecTwo.size() == 3);
SeedingHitSet vecThree(containerThree);
REQUIRE(vecThree.size() == 0);
}

SECTION("Check braced constructor") {
SeedingHitSet listOne({hitThree, hitFive, hitFour});
REQUIRE(listOne.size() == 3);
SeedingHitSet listTwo({HIT_NULL, hitTwo, hitFour});
REQUIRE(listTwo.size() == 0);
SeedingHitSet listThree({hitTwo, hitFour, hitOne, hitFive, HIT_NULL});
REQUIRE(listThree.size() == 4);
}

SECTION("Check two hits constructor") {
SeedingHitSet twoHitsOne(hitOne, hitTwo);
REQUIRE(twoHitsOne.size() == 2);
SeedingHitSet twoHitsTwo(HIT_NULL, hitTwo);
REQUIRE(twoHitsTwo.size() == 0);
SeedingHitSet twoHitsThree(hitOne, HIT_NULL);
REQUIRE(twoHitsThree.size() == 0);
}

SECTION("Check three hits constructor") {
SeedingHitSet threeHitsOne(hitOne, hitTwo, hitThree);
REQUIRE(threeHitsOne.size() == 3);
SeedingHitSet threeHitsTwo(hitOne, hitTwo, HIT_NULL);
REQUIRE(threeHitsTwo.size() == 2);
SeedingHitSet threeHitsThree(hitOne, HIT_NULL, hitThree);
REQUIRE(threeHitsThree.size() == 0);
}
SECTION("Check four hits constructor") {
SeedingHitSet fourHitsOne(hitOne, hitTwo, hitThree, hitFour);
REQUIRE(fourHitsOne.size() == 4);
SeedingHitSet fourHitsTwo(hitOne, hitTwo, HIT_NULL, hitFour);
REQUIRE(fourHitsTwo.size() == 2);
SeedingHitSet fourHitsThree(hitOne, hitTwo, hitThree, HIT_NULL);
REQUIRE(fourHitsThree.size() == 3);
SeedingHitSet fourHitsFour(hitOne, HIT_NULL, hitThree, hitFour);
REQUIRE(fourHitsFour.size() == 0);
}
}
2 changes: 2 additions & 0 deletions RecoTracker/TkSeedingLayers/test/test_catch2_main.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,2 @@
#define CATCH_CONFIG_MAIN // This tells Catch to provide a main() - only do this in one cpp file
#include "catch.hpp"