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

Add pT Runtime Toggle + pT Toggle for Geometry Inputs #39

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
15 changes: 9 additions & 6 deletions RecoTracker/LST/plugins/alpaka/LSTModulesDevESProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -25,29 +25,32 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE {

private:
edm::ESGetToken<SDL::LSTESHostData<SDL::Dev>, TrackerRecoGeometryRecord> lstESHostToken_;
std::string ptCutLabel_;
};

LSTModulesDevESProducer::LSTModulesDevESProducer(const edm::ParameterSet& iConfig) : ESProducer(iConfig) {
setWhatProduced(this, &LSTModulesDevESProducer::produceHost);
auto cc = setWhatProduced(this, &LSTModulesDevESProducer::produceDevice);
lstESHostToken_ = cc.consumes();
LSTModulesDevESProducer::LSTModulesDevESProducer(const edm::ParameterSet& iConfig)
: ESProducer(iConfig), ptCutLabel_(iConfig.getParameter<std::string>("ptCutLabel")) {
setWhatProduced(this, &LSTModulesDevESProducer::produceHost, ptCutLabel_);
auto cc = setWhatProduced(this, &LSTModulesDevESProducer::produceDevice, ptCutLabel_);
lstESHostToken_ = cc.consumes(edm::ESInputTag("", ptCutLabel_));
}

void LSTModulesDevESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<std::string>("ptCutLabel", "0.8");
descriptions.addWithDefaultLabel(desc);
}

std::unique_ptr<SDL::LSTESHostData<SDL::Dev>> LSTModulesDevESProducer::produceHost(
TrackerRecoGeometryRecord const& iRecord) {
return SDL::loadAndFillESHost();
return SDL::loadAndFillESHost(ptCutLabel_);
}

std::unique_ptr<SDL::LSTESDeviceData<SDL::Dev>> LSTModulesDevESProducer::produceDevice(
device::Record<TrackerRecoGeometryRecord> const& iRecord) {
auto const& lstESHostData = iRecord.get(lstESHostToken_);
SDL::QueueAcc& queue = iRecord.queue();
return SDL::loadAndFillESDevice(queue, &lstESHostData);
return SDL::loadAndFillESDevice(queue, &lstESHostData, ptCutLabel_);
}

} // namespace ALPAKA_ACCELERATOR_NAMESPACE
Expand Down
11 changes: 9 additions & 2 deletions RecoTracker/LST/plugins/alpaka/LSTProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -29,8 +29,9 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE {
: lstPixelSeedInputToken_{consumes<LSTPixelSeedInput>(config.getParameter<edm::InputTag>("pixelSeedInput"))},
lstPhase2OTHitsInputToken_{
consumes<LSTPhase2OTHitsInput>(config.getParameter<edm::InputTag>("phase2OTHitsInput"))},
lstESToken_{esConsumes()},
lstESToken_{esConsumes(edm::ESInputTag("", config.getParameter<std::string>("ptCutLabel")))},
verbose_(config.getParameter<bool>("verbose")),
ptCut_(config.getParameter<double>("ptCut")),
nopLSDupClean_(config.getParameter<bool>("nopLSDupClean")),
tcpLSTriplets_(config.getParameter<bool>("tcpLSTriplets")),
lstOutputToken_{produces()} {}
Expand All @@ -44,6 +45,7 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE {

lst_.run(event.queue(),
verbose_,
static_cast<float>(ptCut_),
&lstESDeviceData,
pixelSeeds.px(),
pixelSeeds.py(),
Expand Down Expand Up @@ -81,6 +83,8 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE {
desc.add<edm::InputTag>("pixelSeedInput", edm::InputTag{"lstPixelSeedInputProducer"});
desc.add<edm::InputTag>("phase2OTHitsInput", edm::InputTag{"lstPhase2OTHitsInputProducer"});
desc.add<bool>("verbose", false);
desc.add<double>("ptCut", 0.8);
desc.add<std::string>("ptCutLabel", "0.8");
desc.add<bool>("nopLSDupClean", false);
desc.add<bool>("tcpLSTriplets", false);
descriptions.addWithDefaultLabel(desc);
Expand All @@ -90,7 +94,10 @@ namespace ALPAKA_ACCELERATOR_NAMESPACE {
edm::EDGetTokenT<LSTPixelSeedInput> lstPixelSeedInputToken_;
edm::EDGetTokenT<LSTPhase2OTHitsInput> lstPhase2OTHitsInputToken_;
device::ESGetToken<SDL::LSTESDeviceData<SDL::Dev>, TrackerRecoGeometryRecord> lstESToken_;
const bool verbose_, nopLSDupClean_, tcpLSTriplets_;
const bool verbose_;
const double ptCut_;
const bool nopLSDupClean_;
const bool tcpLSTriplets_;
GNiendorf marked this conversation as resolved.
Show resolved Hide resolved
edm::EDPutTokenT<LSTOutput> lstOutputToken_;

SDL::LST<SDL::Acc> lst_;
Expand Down
2 changes: 1 addition & 1 deletion RecoTracker/LSTCore/BuildFile.xml
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,7 @@
<use name="boost_header"/>
<use name="root"/>
<use name="HeterogeneousCore/AlpakaInterface"/>
<flags CXXFLAGS="-DCACHE_ALLOC -DT4FromT3 -DUSE_RZCHI2 -DUSE_T5_DNN -DPT_CUT=0.8 -DDUP_pLS -DDUP_T5 -DDUP_pT5 -DDUP_pT3 -DCrossclean_T5 -DCrossclean_pT3 -Wshadow"/>
<flags CXXFLAGS="-DCACHE_ALLOC -DT4FromT3 -DUSE_RZCHI2 -DUSE_T5_DNN -DDUP_pLS -DDUP_T5 -DDUP_pT5 -DDUP_pT3 -DCrossclean_T5 -DCrossclean_pT3 -Wshadow"/>
<flags ALPAKA_BACKENDS="1"/>
<export>
<lib name="1"/>
Expand Down
6 changes: 0 additions & 6 deletions RecoTracker/LSTCore/interface/alpaka/Constants.h
Original file line number Diff line number Diff line change
Expand Up @@ -99,11 +99,6 @@ namespace SDL {
return WorkDiv(adjustedBlocks, adjustedThreads, elementsPerThreadArg);
}

// If a compile time flag does not define PT_CUT, default to 0.8 (GeV)
#ifndef PT_CUT
constexpr float PT_CUT = 0.8f;
#endif

const unsigned int MAX_BLOCKS = 80;
const unsigned int MAX_CONNECTED_MODULES = 40;

Expand All @@ -129,7 +124,6 @@ namespace SDL {
ALPAKA_STATIC_ACC_MEM_GLOBAL const float k2Rinv1GeVf = (2.99792458e-3 * 3.8) / 2;
ALPAKA_STATIC_ACC_MEM_GLOBAL const float kR1GeVf = 1. / (2.99792458e-3 * 3.8);
ALPAKA_STATIC_ACC_MEM_GLOBAL const float sinAlphaMax = 0.95;
ALPAKA_STATIC_ACC_MEM_GLOBAL const float ptCut = PT_CUT;
ALPAKA_STATIC_ACC_MEM_GLOBAL const float deltaZLum = 15.0;
ALPAKA_STATIC_ACC_MEM_GLOBAL const float pixelPSZpitch = 0.15;
ALPAKA_STATIC_ACC_MEM_GLOBAL const float strip2SZpitch = 5.0;
Expand Down
4 changes: 3 additions & 1 deletion RecoTracker/LSTCore/interface/alpaka/LST.h
Original file line number Diff line number Diff line change
Expand Up @@ -23,6 +23,7 @@ namespace SDL {

void run(QueueAcc& queue,
bool verbose,
const float ptCut,
const LSTESDeviceData<Dev>* deviceESData,
const std::vector<float> see_px,
const std::vector<float> see_py,
Expand Down Expand Up @@ -69,7 +70,8 @@ namespace SDL {
const std::vector<unsigned int> ph2_detId,
const std::vector<float> ph2_x,
const std::vector<float> ph2_y,
const std::vector<float> ph2_z);
const std::vector<float> ph2_z,
const float ptCut);

void getOutput(SDL::Event<Acc>& event);
std::vector<unsigned int> getHitIdxs(const short trackCandidateType,
Expand Down
6 changes: 4 additions & 2 deletions RecoTracker/LSTCore/interface/alpaka/LSTESData.h
Original file line number Diff line number Diff line change
Expand Up @@ -72,8 +72,10 @@ namespace SDL {
pixelMapping(pixelMappingIn) {}
};

std::unique_ptr<LSTESHostData<Dev>> loadAndFillESHost();
std::unique_ptr<LSTESDeviceData<Dev>> loadAndFillESDevice(SDL::QueueAcc& queue, const LSTESHostData<Dev>* hostData);
std::unique_ptr<LSTESHostData<Dev>> loadAndFillESHost(std::string& ptCutLabel);
std::unique_ptr<LSTESDeviceData<Dev>> loadAndFillESDevice(SDL::QueueAcc& queue,
const LSTESHostData<Dev>* hostData,
std::string& ptCutLabel);

} // namespace SDL

Expand Down
18 changes: 12 additions & 6 deletions RecoTracker/LSTCore/src/alpaka/Event.dev.cc
Original file line number Diff line number Diff line change
Expand Up @@ -440,7 +440,8 @@ void SDL::Event<SDL::Acc>::createMiniDoublets() {
*modulesBuffers_->data(),
*hitsInGPU,
*mdsInGPU,
*rangesInGPU));
*rangesInGPU,
ptCut));

alpaka::enqueue(queue, createMiniDoubletsInGPUv2Task);

Expand Down Expand Up @@ -485,7 +486,8 @@ void SDL::Event<SDL::Acc>::createSegmentsWithModuleMap() {
*modulesBuffers_->data(),
*mdsInGPU,
*segmentsInGPU,
*rangesInGPU));
*rangesInGPU,
ptCut));

alpaka::enqueue(queue, createSegmentsInGPUv2Task);

Expand Down Expand Up @@ -593,7 +595,8 @@ void SDL::Event<SDL::Acc>::createTriplets() {
*tripletsInGPU,
*rangesInGPU,
alpaka::getPtrNative(index_gpu_buf),
nonZeroModules));
nonZeroModules,
ptCut));

alpaka::enqueue(queue, createTripletsInGPUv2Task);

Expand Down Expand Up @@ -872,7 +875,8 @@ void SDL::Event<SDL::Acc>::createPixelTriplets() {
*pixelTripletsInGPU,
alpaka::getPtrNative(connectedPixelSize_dev_buf),
alpaka::getPtrNative(connectedPixelIndex_dev_buf),
nInnerSegments));
nInnerSegments,
ptCut));

alpaka::enqueue(queue, createPixelTripletsInGPUFromMapv2Task);
alpaka::wait(queue);
Expand Down Expand Up @@ -951,7 +955,8 @@ void SDL::Event<SDL::Acc>::createQuintuplets() {
*tripletsInGPU,
*quintupletsInGPU,
*rangesInGPU,
nEligibleT5Modules));
nEligibleT5Modules,
ptCut));

alpaka::enqueue(queue, createQuintupletsInGPUv2Task);

Expand Down Expand Up @@ -1102,7 +1107,8 @@ void SDL::Event<SDL::Acc>::createPixelQuintuplets() {
alpaka::getPtrNative(connectedPixelSize_dev_buf),
alpaka::getPtrNative(connectedPixelIndex_dev_buf),
nInnerSegments,
*rangesInGPU));
*rangesInGPU,
ptCut));

alpaka::enqueue(queue, createPixelQuintupletsInGPUFromMapv2Task);

Expand Down
8 changes: 7 additions & 1 deletion RecoTracker/LSTCore/src/alpaka/Event.h
Original file line number Diff line number Diff line change
Expand Up @@ -25,6 +25,7 @@ namespace SDL {
class Event<SDL::Acc> {
private:
QueueAcc queue;
const float ptCut;
Dev devAcc;
DevHost devHost;
bool addObjects;
Expand Down Expand Up @@ -91,8 +92,9 @@ namespace SDL {
public:
// Constructor used for CMSSW integration. Uses an external queue.
template <typename TQueue>
Event(bool verbose, TQueue const& q, const LSTESDeviceData<Dev>* deviceESData)
Event(bool verbose, const float pt_cut, TQueue const& q, const LSTESDeviceData<Dev>* deviceESData)
: queue(q),
ptCut(pt_cut),
devAcc(alpaka::getDev(q)),
devHost(cms::alpakatools::host()),
nModules_(deviceESData->nModules),
Expand All @@ -101,6 +103,10 @@ namespace SDL {
modulesBuffers_(deviceESData->modulesBuffers),
pixelMapping_(deviceESData->pixelMapping),
endcapGeometry_(deviceESData->endcapGeometry) {
if (pt_cut < 0.6f) {
throw std::invalid_argument("Minimum pT cut must be at least 0.6 GeV. Provided value: " +
std::to_string(pt_cut));
}
init(verbose);
}
void resetEvent();
Expand Down
13 changes: 8 additions & 5 deletions RecoTracker/LSTCore/src/alpaka/LST.dev.cc
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@ using XYZVector = ROOT::Math::XYZVector;

void SDL::LST<SDL::Acc>::run(SDL::QueueAcc& queue,
bool verbose,
const float ptCut,
const LSTESDeviceData<SDL::Dev>* deviceESData,
const std::vector<float> see_px,
const std::vector<float> see_py,
Expand All @@ -29,7 +30,7 @@ void SDL::LST<SDL::Acc>::run(SDL::QueueAcc& queue,
const std::vector<float> ph2_z,
bool no_pls_dupclean,
bool tc_pls_triplets) {
auto event = SDL::Event<Acc>(verbose, queue, deviceESData);
auto event = SDL::Event<Acc>(verbose, ptCut, queue, deviceESData);
prepareInput(see_px,
see_py,
see_pz,
Expand All @@ -48,7 +49,8 @@ void SDL::LST<SDL::Acc>::run(SDL::QueueAcc& queue,
ph2_detId,
ph2_x,
ph2_y,
ph2_z);
ph2_z,
ptCut);

event.addHitToEvent(in_trkX_, in_trkY_, in_trkZ_, in_hitId_, in_hitIdxs_);
event.addPixelSegmentToEvent(in_hitIndices_vec0_,
Expand Down Expand Up @@ -188,7 +190,8 @@ void SDL::LST<SDL::Acc>::prepareInput(const std::vector<float> see_px,
const std::vector<unsigned int> ph2_detId,
const std::vector<float> ph2_x,
const std::vector<float> ph2_y,
const std::vector<float> ph2_z) {
const std::vector<float> ph2_z,
const float ptCut) {
unsigned int count = 0;
auto n_see = see_stateTrajGlbPx.size();
std::vector<float> px_vec;
Expand Down Expand Up @@ -240,7 +243,7 @@ void SDL::LST<SDL::Acc>::prepareInput(const std::vector<float> see_px,
float eta = p3LH.eta();
float ptErr = see_ptErr[iSeed];

if ((ptIn > 0.8 - 2 * ptErr)) {
if ((ptIn > ptCut - 2 * ptErr)) {
XYZVector r3LH(see_stateTrajGlbX[iSeed], see_stateTrajGlbY[iSeed], see_stateTrajGlbZ[iSeed]);
XYZVector p3PCA(see_px[iSeed], see_py[iSeed], see_pz[iSeed]);
XYZVector r3PCA(calculateR3FromPCA(p3PCA, see_dxy[iSeed], see_dz[iSeed]));
Expand All @@ -256,7 +259,7 @@ void SDL::LST<SDL::Acc>::prepareInput(const std::vector<float> see_px,

if (ptIn >= 2.0)
pixtype = 0;
else if (ptIn >= (0.8 - 2 * ptErr) and ptIn < 2.0) {
else if (ptIn >= (ptCut - 2 * ptErr) and ptIn < 2.0) {
if (pixelSegmentDeltaPhiChange >= 0)
pixtype = 1;
else
Expand Down
28 changes: 15 additions & 13 deletions RecoTracker/LSTCore/src/alpaka/LSTESData.dev.cc
Original file line number Diff line number Diff line change
Expand Up @@ -41,21 +41,22 @@ namespace {
void loadMapsHost(SDL::MapPLStoLayer& pLStoLayer,
std::shared_ptr<SDL::EndcapGeometryHost<SDL::Dev>> endcapGeometry,
std::shared_ptr<SDL::TiltedGeometry<SDL::Dev>> tiltedGeometry,
std::shared_ptr<SDL::ModuleConnectionMap<SDL::Dev>> moduleConnectionMap) {
std::shared_ptr<SDL::ModuleConnectionMap<SDL::Dev>> moduleConnectionMap,
std::string& ptCutLabel) {
// Module orientation information (DrDz or phi angles)
auto endcap_geom =
get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt0.8/endcap_orientation.bin");
auto tilted_geom = get_absolute_path_after_check_file_exists(
trackLooperDir() + "/data/OT800_IT615_pt0.8/tilted_barrel_orientation.bin");
auto endcap_geom = get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt" +
ptCutLabel + "/endcap_orientation.bin");
auto tilted_geom = get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt" +
ptCutLabel + "/tilted_barrel_orientation.bin");
// Module connection map (for line segment building)
auto mappath = get_absolute_path_after_check_file_exists(
trackLooperDir() + "/data/OT800_IT615_pt0.8/module_connection_tracing_merged.bin");
auto mappath = get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt" + ptCutLabel +
"/module_connection_tracing_merged.bin");

endcapGeometry->load(endcap_geom);
tiltedGeometry->load(tilted_geom);
moduleConnectionMap->load(mappath);

auto pLSMapDir = trackLooperDir() + "/data/OT800_IT615_pt0.8/pixelmap/pLS_map";
auto pLSMapDir = trackLooperDir() + "/data/OT800_IT615_pt" + ptCutLabel + "/pixelmap/pLS_map";
const std::array<std::string, 4> connects{
{"_layer1_subdet5", "_layer2_subdet5", "_layer1_subdet4", "_layer2_subdet4"}};
std::string path;
Expand All @@ -76,17 +77,18 @@ namespace {
}
} // namespace

std::unique_ptr<SDL::LSTESHostData<SDL::Dev>> SDL::loadAndFillESHost() {
std::unique_ptr<SDL::LSTESHostData<SDL::Dev>> SDL::loadAndFillESHost(std::string& ptCutLabel) {
auto pLStoLayer = std::make_shared<SDL::MapPLStoLayer>();
auto endcapGeometry = std::make_shared<SDL::EndcapGeometryHost<SDL::Dev>>();
auto tiltedGeometry = std::make_shared<SDL::TiltedGeometry<SDL::Dev>>();
auto moduleConnectionMap = std::make_shared<SDL::ModuleConnectionMap<SDL::Dev>>();
::loadMapsHost(*pLStoLayer, endcapGeometry, tiltedGeometry, moduleConnectionMap);
::loadMapsHost(*pLStoLayer, endcapGeometry, tiltedGeometry, moduleConnectionMap, ptCutLabel);
return std::make_unique<LSTESHostData<SDL::Dev>>(pLStoLayer, endcapGeometry, tiltedGeometry, moduleConnectionMap);
}

std::unique_ptr<SDL::LSTESDeviceData<SDL::Dev>> SDL::loadAndFillESDevice(SDL::QueueAcc& queue,
const LSTESHostData<SDL::Dev>* hostData) {
const LSTESHostData<SDL::Dev>* hostData,
std::string& ptCutLabel) {
SDL::Dev const& devAccIn = alpaka::getDev(queue);
uint16_t nModules;
uint16_t nLowerModules;
Expand All @@ -96,8 +98,8 @@ std::unique_ptr<SDL::LSTESDeviceData<SDL::Dev>> SDL::loadAndFillESDevice(SDL::Qu
auto pixelMapping = std::make_shared<SDL::pixelMap>();
auto moduleConnectionMap = hostData->moduleConnectionMap;

auto path =
get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt0.8/sensor_centroids.bin");
auto path = get_absolute_path_after_check_file_exists(trackLooperDir() + "/data/OT800_IT615_pt" + ptCutLabel +
"/sensor_centroids.bin");
SDL::loadModulesFromFile(queue,
hostData->mapPLStoLayer.get(),
path.c_str(),
Expand Down
Loading