ATLAS Offline Software
Loading...
Searching...
No Matches
columnar::ColumnarPhysLiteTest Struct Reference

#include <ColumnarPhysliteTest.h>

Inheritance diagram for columnar::ColumnarPhysLiteTest:
Collaboration diagram for columnar::ColumnarPhysLiteTest:

Public Member Functions

 ColumnarPhysLiteTest ()
 ~ColumnarPhysLiteTest ()
std::string makeUniqueName ()
 make a unique tool name to be used in unit tests
void setupKnownColumns (std::span< const TestUtils::TestDefinition > testDefinitions)
void setupColumns (const ColumnVectorHeader &columnHeader)
void doCall (const TestUtils::TestDefinition &testDefinition)
void doCallMulti (const std::vector< TestUtils::TestDefinition > &testDefinitions)

Static Public Member Functions

static bool checkMode ()
 check whether we have the right mode

Public Attributes

std::unique_ptr< TFile > file
TTree * tree = nullptr
std::unique_ptr< ROOT::RNTupleReader > rntreader
std::unique_ptr< ROOT::Experimental::RNTupleInspector > inspector
columnar::TestUtils::RNTupleBackendrntbackend = nullptr
std::vector< std::shared_ptr< TestUtils::IColumnData > > knownColumns
std::vector< std::shared_ptr< TestUtils::IColumnData > > usedColumns
std::unordered_map< std::string, const std::vector< ColumnarOffsetType > * > offsetColumns

Detailed Description

Definition at line 42 of file ColumnarPhysliteTest.h.

Constructor & Destructor Documentation

◆ ColumnarPhysLiteTest()

columnar::ColumnarPhysLiteTest::ColumnarPhysLiteTest ( )

Definition at line 2141 of file ColumnarPhysliteTest.cxx.

2143 {
2144 static std::once_flag flag;
2145 std::call_once (flag, [] ()
2146 {
2147#ifdef XAOD_STANDALONE
2148 xAOD::Init().ignore();
2149#else
2150 POOL::Init();
2151#endif
2152 });
2153
2154 auto userConfiguration = TestUtils::UserConfiguration::fromEnvironment();
2155 if (userConfiguration.isrntuple)
2156 {
2157 auto* fileName = getenv("ASG_TEST_FILE_RNTUPLE_LITE_MC");
2158 if (fileName == nullptr)
2159 throw std::runtime_error("missing ASG_TEST_FILE_RNTUPLE_LITE_MC");
2160 rntreader = ROOT::RNTupleReader::Open("EventData", fileName);
2161 inspector = ROOT::Experimental::RNTupleInspector::Create("EventData", fileName);
2162 rntbackend = new TestUtils::RNTupleBackend{rntreader.get(), inspector.get()};
2163 if (!rntreader or !inspector)
2164 throw std::runtime_error("failed to open rntuple");
2165 } else
2166 {
2167 auto* fileName = getenv("ASG_TEST_FILE_LITE_MC");
2168 if (fileName == nullptr)
2169 throw std::runtime_error("missing ASG_TEST_FILE_LITE_MC");
2170 file.reset(TFile::Open(fileName, "READ"));
2171 if (!file)
2172 throw std::runtime_error("failed to open file");
2173 tree = dynamic_cast<TTree*>(file->Get("CollectionTree"));
2174 if (!tree)
2175 throw std::runtime_error("failed to open rntuple");
2176 }
2177 }
IAppMgrUI * Init(const char *options="POOLRootAccess/basic.opts")
Bootstraps (creates and configures) the Gaudi Application with the provided options file.
std::string getenv(const std::string &variableName)
get an environment variable
bool flag
Definition master.py:29
StatusCode Init(const char *appname)
Function initialising ROOT/PyROOT for using the ATLAS EDM.
Definition Init.cxx:31
columnar::TestUtils::RNTupleBackend * rntbackend
std::unique_ptr< ROOT::Experimental::RNTupleInspector > inspector
std::unique_ptr< ROOT::RNTupleReader > rntreader
static UserConfiguration fromEnvironment()
create a UserConfiguration, loading from the file pointed to by the COLUMNAR_TEST_CONFIG environment ...

◆ ~ColumnarPhysLiteTest()

columnar::ColumnarPhysLiteTest::~ColumnarPhysLiteTest ( )

Definition at line 2179 of file ColumnarPhysliteTest.cxx.

2180 {
2181 if (rntbackend)
2182 delete rntbackend;
2183 }

Member Function Documentation

◆ checkMode()

bool columnar::ColumnarPhysLiteTest::checkMode ( )
static

check whether we have the right mode

Definition at line 2191 of file ColumnarPhysliteTest.cxx.

2193 {
2194 return true;
2195 }

◆ doCall()

void columnar::ColumnarPhysLiteTest::doCall ( const TestUtils::TestDefinition & testDefinition)

Definition at line 2625 of file ColumnarPhysliteTest.cxx.

2626 {
2627 doCallMulti ({testDefinition});
2628 }
void doCallMulti(const std::vector< TestUtils::TestDefinition > &testDefinitions)

◆ doCallMulti()

void columnar::ColumnarPhysLiteTest::doCallMulti ( const std::vector< TestUtils::TestDefinition > & testDefinitions)

Definition at line 2631 of file ColumnarPhysliteTest.cxx.

2632 {
2633 using namespace asg::msgUserCode;
2634 auto userConfiguration = TestUtils::UserConfiguration::fromEnvironment();
2635
2636 // apply systematics for all test definitions
2637 for (const auto& td : testDefinitions) {
2638 if (!td.sysName.empty()) {
2639 auto* sysTool = dynamic_cast<CP::ISystematicsTool*>(td.tool);
2640 if (!sysTool)
2641 throw std::runtime_error("tool does not support systematics");
2642 std::cout << "applying systematic variation: " << td.sysName << std::endl;
2643 if (sysTool->applySystematicVariation(CP::SystematicSet(td.sysName))
2644 .isFailure())
2645 throw std::runtime_error("failed to apply systematic variation: " +
2646 td.sysName);
2647 }
2648 }
2649
2650 if constexpr (columnarAccessMode == 2) {
2651 // Create shared column header for all tools
2652 ColumnVectorHeader columnHeader;
2653
2654 // Build vector of ToolData from all testDefinitions
2655 std::vector<TestUtils::ToolData> toolDataVec;
2656 for (const auto& td : testDefinitions)
2657 toolDataVec.emplace_back(userConfiguration, td, columnHeader);
2658
2659 setupKnownColumns(testDefinitions);
2660 // Set up columns using the shared header (all tools have already
2661 // registered their columns via ToolColumnVectorMap, so we get all columns
2662 // from the header)
2663 setupColumns(columnHeader);
2664
2665 // connect column indices from header to each column for direct setting
2666 for (auto& column : usedColumns)
2667 column->connectColumnIndices(columnHeader);
2668
2669 Benchmark benchmarkEmpty("empty");
2670 Benchmark benchmarkCheck("", userConfiguration.batchSize);
2671 auto numberOfEvents = 0;
2672 if (tree) {
2673 numberOfEvents = tree->GetEntries();
2674 } else if (rntbackend) {
2675 numberOfEvents = rntreader->GetNEntries();
2676 }
2677 Long64_t entry = 0;
2678 const auto startTime = std::chrono::high_resolution_clock::now();
2679 bool endLoop = false;
2680 for (; !endLoop; ++entry) {
2681 // just sample how much overhead there is for starting and
2682 // stopping the timer
2683 benchmarkEmpty.startTimer();
2684 benchmarkEmpty.stopTimer();
2685 ColumnVectorData columnData(&columnHeader);
2686 for (auto& column : usedColumns)
2687 column->getEntry(entry % numberOfEvents);
2688 if ((entry + 1) % userConfiguration.batchSize == 0) {
2689 if (entry < numberOfEvents) {
2690 for (auto& column : usedColumns)
2691 column->collectColumnData();
2692 }
2693 for (auto& column : usedColumns)
2694 column->setData(columnData);
2695
2696 // Check data once (shared column data)
2697 benchmarkCheck.startTimer();
2698 columnData.checkData();
2699 benchmarkCheck.stopTimer();
2700 // Call each tool
2701 for (auto& toolData : toolDataVec) {
2702 toolData.call(columnData);
2703 }
2704 for (auto& column : usedColumns)
2705 column->clearColumns();
2706 if ((std::chrono::high_resolution_clock::now() - startTime) >
2707 userConfiguration.targetTime)
2708 endLoop = true;
2709 } else if (entry + 1 == numberOfEvents) {
2710 for (auto& column : usedColumns)
2711 column->collectColumnData();
2712 }
2713 }
2714 std::cout << "Entries in file: " << numberOfEvents << std::endl;
2715 std::cout << "Total entries read: " << entry << std::endl;
2716 const float emptyTime = benchmarkEmpty.getEntryTime(0).value();
2717 std::cout << "Empty benchmark time: " << emptyTime << "ns (tick=" << Benchmark::getTickDuration() << "ns)" << std::endl;
2718 benchmarkEmpty.setSilence();
2719 const auto checkTime = benchmarkCheck.getEntryTime(emptyTime);
2720 if (checkTime)
2721 std::cout << "Check data time: " << checkTime.value() << "ns" << std::endl;
2722 benchmarkCheck.setSilence();
2723 {
2724 std::vector<TestUtils::BranchPerfData> branchPerfData;
2725 TestUtils::BranchPerfData summary;
2726 summary.name = "total";
2727 summary.timeRead = 0;
2728 summary.timeUnpack = 0;
2729 summary.timeShallowCopy = 0;
2730 summary.entrySize = 0;
2731 summary.uncompressedSize = 0;
2732 summary.numBaskets = 0;
2733 summary.entries = std::nullopt;
2734 summary.nullEntries = std::nullopt;
2735 for (auto& column : usedColumns)
2736 {
2737 branchPerfData.push_back (column->getPerfData (emptyTime));
2738 summary.timeRead.value() += branchPerfData.back().timeRead.value_or(0);
2739 summary.timeUnpack.value() += branchPerfData.back().timeUnpack.value_or(0);
2740 summary.entrySize.value() += branchPerfData.back().entrySize.value_or(0);
2741 summary.uncompressedSize.value() += branchPerfData.back().uncompressedSize.value_or(0);
2742 summary.numBaskets.value() += branchPerfData.back().numBaskets.value_or(0);
2743 summary.timeShallowCopy.value() += branchPerfData.back().timeShallowCopy.value_or(0);
2744 }
2745 std::sort (branchPerfData.begin(), branchPerfData.end(), [] (const auto& a, const auto& b) {return a.name < b.name;});
2746 branchPerfData.insert (branchPerfData.end(), summary);
2747 const std::size_t nameWidth = std::max_element (branchPerfData.begin(), branchPerfData.end(), [] (const auto& a, const auto& b) {return a.name.size() < b.name.size();})->name.size();
2748 std::string label = userConfiguration.isrntuple ? "field name" : "branch name";
2749 std::string header = std::format ("{:{}} | read(ns) | unpack(ns) | size(B) | rate(MB/s) | compression | baskets | entries | null", label, nameWidth);
2750 std::cout << "\n" << header << std::endl;
2751 std::cout << std::string (header.size(), '-') << std::endl;
2752 for (auto& data : branchPerfData)
2753 {
2754 if (data.name == "total")
2755 std::cout << std::string (header.size(), '-') << std::endl;
2756 std::cout << std::format ("{:{}} |", data.name, nameWidth);
2757 if (data.timeRead)
2758 std::cout << std::format ("{:>9.0f} |", data.timeRead.value());
2759 else
2760 std::cout << " |";
2761 if (data.timeUnpack)
2762 std::cout << std::format ("{:>11.1f} |", data.timeUnpack.value());
2763 else
2764 std::cout << " |";
2765 if (data.entrySize)
2766 std::cout << std::format ("{:>8.1f} |", data.entrySize.value());
2767 else
2768 std::cout << " |";
2769 if (data.timeRead && data.entrySize)
2770 std::cout << std::format ("{:>11.1f} |", (data.entrySize.value() / (data.timeRead.value() * 1e-3 * 1.024 * 1.024)));
2771 else
2772 std::cout << " |";
2773 if (data.entrySize && data.uncompressedSize)
2774 std::cout << std::format ("{:>12.2f} |", float (data.uncompressedSize.value()) / data.entrySize.value());
2775 else
2776 std::cout << " |";
2777 if (data.numBaskets)
2778 std::cout << std::format ("{:>8} |", data.numBaskets.value());
2779 else
2780 std::cout << " |";
2781 if (data.entries)
2782 std::cout << std::format ("{:>8.2f} |", static_cast<float>(data.entries.value())/numberOfEvents);
2783 else
2784 std::cout << " |";
2785 if (data.nullEntries && data.entries)
2786 std::cout << std::format ("{:>4.0f}%", static_cast<float>(data.nullEntries.value()) / data.entries.value() * 100.0f);
2787 std::cout << std::endl;
2788 }
2789 }
2790 {
2791 std::vector<TestUtils::ToolPerfData> toolPerfData;
2792 for (auto& toolData : toolDataVec)
2793 {
2794 toolPerfData.emplace_back ();
2795 toolPerfData.back().name = toolData.name;
2796 toolPerfData.back().timeCall = toolData.benchmarkCall.getEntryTime (emptyTime);
2797 if (userConfiguration.runToolTwice)
2798 toolPerfData.back().timeCall2 = toolData.benchmarkCall2.getEntryTime (emptyTime);
2799 }
2800 const std::size_t nameWidth = std::max_element (toolPerfData.begin(), toolPerfData.end(), [] (const auto& a, const auto& b) {return a.name.size() < b.name.size();})->name.size();
2801 std::string header = std::format ("{:{}} | call(ns) | call2(ns)", "tool name", nameWidth);
2802 std::cout << "\n" << header << std::endl;
2803 std::cout << std::string (header.size(), '-') << std::endl;
2804 for (auto& data : toolPerfData)
2805 {
2806 std::cout << std::format ("{:{}} |", data.name, nameWidth);
2807 if (data.timeCall)
2808 std::cout << std::format ("{:>9.0f} |", data.timeCall.value());
2809 else
2810 std::cout << " |";
2811 if (data.timeCall2)
2812 std::cout << std::format ("{:>10.0f}", data.timeCall2.value());
2813 else
2814 std::cout << " ";
2815 std::cout << std::endl;
2816 }
2817 // Add totals line for multiple tools
2818 if (toolPerfData.size() > 1)
2819 {
2820 std::optional<float> totalCall, totalCall2;
2821 for (const auto& data : toolPerfData)
2822 {
2823 if (data.timeCall)
2824 totalCall = totalCall.value_or (0) + data.timeCall.value();
2825 if (data.timeCall2)
2826 totalCall2 = totalCall2.value_or (0) + data.timeCall2.value();
2827 }
2828 std::cout << std::string (header.size(), '-') << std::endl;
2829 std::cout << std::format ("{:{}} |", "total", nameWidth);
2830 if (totalCall)
2831 std::cout << std::format ("{:>9.0f} |", totalCall.value());
2832 else
2833 std::cout << " |";
2834 if (totalCall2)
2835 std::cout << std::format ("{:>10.0f}", totalCall2.value());
2836 else
2837 std::cout << " ";
2838 std::cout << std::endl;
2839 }
2840 }
2841 } else if constexpr (columnarAccessMode == 0)
2842 {
2843 TestUtils::runXaodTest (userConfiguration, testDefinitions, file.get());
2844 } else if constexpr (columnarAccessMode == 100)
2845 {
2846 const auto& testDefinition = testDefinitions[0];
2847 TestUtils::runXaodArrayTest (userConfiguration, testDefinition, file.get());
2848 }
2849 }
void checkTime()
int numberOfEvents()
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
static Double_t a
static float getTickDuration()
Definition Benchmark.h:86
std::string label(const std::string &format, int i)
Definition label.h:19
void runXaodArrayTest(const UserConfiguration &userConfiguration, const TestDefinition &testDefinition, TFile *file)
void runXaodTest(const UserConfiguration &userConfiguration, std::span< const TestDefinition > testDefinitions, TFile *file)
constexpr unsigned columnarAccessMode
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
std::vector< std::shared_ptr< TestUtils::IColumnData > > usedColumns
void setupKnownColumns(std::span< const TestUtils::TestDefinition > testDefinitions)
void setupColumns(const ColumnVectorHeader &columnHeader)

◆ makeUniqueName()

std::string columnar::ColumnarPhysLiteTest::makeUniqueName ( )

make a unique tool name to be used in unit tests

Definition at line 2185 of file ColumnarPhysliteTest.cxx.

2186 {
2187 static std::atomic<unsigned> index = 0;
2188 return "UniquePhysliteTestTool" + std::to_string(++index);
2189 }
str index
Definition DeMoScan.py:362

◆ setupColumns()

void columnar::ColumnarPhysLiteTest::setupColumns ( const ColumnVectorHeader & columnHeader)

Definition at line 2517 of file ColumnarPhysliteTest.cxx.

2518 {
2519 using namespace asg::msgUserCode;
2520
2521 // Get all column info directly from the header (all tools have already
2522 // registered their columns via ToolColumnVectorMap)
2523 auto requestedColumns = columnHeader.getAllColumnInfo();
2524
2525 // Print requested columns
2526 for (auto& [name, info] : requestedColumns)
2527 std::cout << "requested columns: " << name << std::endl;
2528
2529 for (auto& column : knownColumns)
2530 {
2531 if (tree)
2532 {
2533 if (column->connect (tree, offsetColumns, requestedColumns))
2534 usedColumns.push_back (column);
2535 } else if (rntbackend)
2536 {
2537 if (column->connect(rntbackend, offsetColumns, requestedColumns))
2538 usedColumns.push_back(column);
2539 }
2540 }
2541
2542 std::set<std::string> unclaimedColumns;
2543 for (auto& column : requestedColumns)
2544 {
2545 if (!column.second.isOptional)
2546 unclaimedColumns.insert (column.first);
2547 else
2548 std::cout << "optional column not claimed: " << column.first << std::endl;
2549 }
2550 std::erase_if (unclaimedColumns, [&] (auto& columnName)
2551 {
2552 const auto& info = requestedColumns.at (columnName);
2553 if (info.accessMode != ColumnAccessMode::output || !info.fixedDimensions.empty())
2554 return false;
2555 auto offsetIter = std::find_if (usedColumns.begin(), usedColumns.end(), [&] (const std::shared_ptr<TestUtils::IColumnData>& column)
2556 {
2557 for (auto& output : column->outputColumns)
2558 {
2559 if (output.name == info.offsetName)
2560 return true;
2561 }
2562 return false;
2563 });
2564 if (offsetIter == usedColumns.end())
2565 return false;
2566 std::shared_ptr<TestUtils::IColumnData> myColumn;
2567 if (tree)
2568 {
2569 if (*info.type == typeid(float))
2570 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<float, BranchReader>>(info.name, 0);
2571 else if (*info.type == typeid(char))
2572 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<char, BranchReader>>(info.name, 0);
2573 else if (*info.type == typeid(std::uint16_t))
2574 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<std::uint16_t, BranchReader>>(info.name, 0);
2575 else if (*info.type == typeid(std::uint64_t))
2576 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<std::uint64_t, BranchReader>>(info.name, 0);
2577 else
2578 {
2579 ANA_MSG_WARNING("unhandled column type: " << info.name << " "<< info.type->name());
2580 return false;
2581 }
2582 } else if (rntbackend)
2583 {
2584 if (*info.type == typeid(float))
2585 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<float, RNTFieldReader>>(info.name,0);
2586 else if (*info.type == typeid(char))
2587 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<char, RNTFieldReader>>(info.name, 0);
2588 else if (*info.type == typeid(std::uint16_t))
2589 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<std::uint16_t, RNTFieldReader>>(info.name, 0);
2590 else if (*info.type == typeid(std::uint64_t))
2591 myColumn = std::make_shared<TestUtils::ColumnDataOutVector<std::uint64_t, RNTFieldReader>>(info.name, 0);
2592 else
2593 {
2594 ANA_MSG_WARNING("unhandled column type: " << info.name << " " << info.type->name());
2595 return false;
2596 }
2597 }
2598 knownColumns.push_back(myColumn);
2599 if (tree) {
2600 if (!myColumn->connect(tree, offsetColumns, requestedColumns))
2601 {
2602 ANA_MSG_WARNING("failed to connect dynamic output column: " << info.name);
2603 return false;
2604 }
2605 } else if (rntbackend)
2606 {
2607 if (!myColumn->connect(rntbackend, offsetColumns, requestedColumns))
2608 {
2609 ANA_MSG_WARNING("failed to connect dynamic output column: " << info.name);
2610 return false;
2611 }
2612 }
2613 usedColumns.push_back(myColumn);
2614 return true;
2615 });
2616 if (!unclaimedColumns.empty())
2617 {
2618 std::string message = "columns not claimed:";
2619 for (auto& column : unclaimedColumns)
2620 message += " " + column;
2621 throw std::runtime_error(message);
2622 }
2623 }
#define ANA_MSG_WARNING(xmsg)
Macro printing warning messages.
@ output
an output column
Definition ColumnInfo.h:24
std::size_t erase_if(T_container &container, T_Func pred)
std::vector< std::shared_ptr< TestUtils::IColumnData > > knownColumns
std::unordered_map< std::string, const std::vector< ColumnarOffsetType > * > offsetColumns

◆ setupKnownColumns()

void columnar::ColumnarPhysLiteTest::setupKnownColumns ( std::span< const TestUtils::TestDefinition > testDefinitions)

Definition at line 2197 of file ColumnarPhysliteTest.cxx.

2198 {
2199 using namespace TestUtils;
2200
2201 knownColumns.push_back (std::make_shared<ColumnDataEventCount> ());
2202
2203 if (tree)
2204 {
2205 tree->SetMakeClass(1);
2206 {
2207 std::unordered_map<std::string, TBranch*> branches;
2208 {
2209 TIter branchIter(tree->GetListOfBranches());
2210 TObject* obj = nullptr;
2211 while ((obj = branchIter()))
2212 {
2213 TBranch* branch = nullptr;
2214 if ((branch = dynamic_cast<TBranch*>(obj)))
2215 {
2216 branches.emplace(branch->GetName(), branch);
2217 TIter subBranchIter(branch->GetListOfBranches());
2218 while ((obj = subBranchIter()))
2219 {
2220 if (auto subBranch = dynamic_cast<TBranch*>(obj))
2221 branches.emplace(subBranch->GetName(), subBranch);
2222 }
2223 }
2224 }
2225 }
2226
2227 for (const auto& [name, branch] : branches)
2228 {
2229 if (name.find("AuxDyn.") != std::string::npos ||
2230 name.find("Aux.") != std::string::npos)
2231 {
2232 TClass* branchClass = nullptr;
2233 EDataType branchType{};
2234 branch->GetExpectedType(branchClass, branchType);
2235 if (branchClass == nullptr)
2236 {
2237 switch (branchType)
2238 {
2239 case kInt_t:
2240 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::int32_t, BranchReader>>(branch->GetName()));
2241 break;
2242 case kUInt_t:
2243 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::uint32_t, BranchReader>>(branch->GetName()));
2244 break;
2245 case kULong_t:
2246 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::uint64_t, BranchReader>>(branch->GetName()));
2247 break;
2248 case kULong64_t:
2249 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::uint64_t, BranchReader>>(branch->GetName()));
2250 break;
2251 case kFloat_t:
2252 knownColumns.push_back(std::make_shared<ColumnDataScalar<float, BranchReader>>(branch->GetName()));
2253 break;
2254 default:
2255 // no-op
2256 break;
2257 }
2258 } else
2259 {
2260 if (*branchClass->GetTypeInfo() == typeid(std::vector<float>))
2261 {
2262 knownColumns.push_back (std::make_shared<ColumnDataVector<float,BranchReader>> (branch->GetName()));
2263 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<char>))
2264 {
2265 knownColumns.push_back (std::make_shared<ColumnDataVector<char,BranchReader>> (branch->GetName()));
2266 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::int8_t>))
2267 {
2268 knownColumns.push_back (std::make_shared<ColumnDataVector<std::int8_t,BranchReader>> (branch->GetName()));
2269 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::uint8_t>))
2270 {
2271 knownColumns.push_back (std::make_shared<ColumnDataVector<std::uint8_t,BranchReader>> (branch->GetName()));
2272 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::int16_t>))
2273 {
2274 knownColumns.push_back (std::make_shared<ColumnDataVector<std::int16_t,BranchReader>> (branch->GetName()));
2275 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::uint16_t>))
2276 {
2277 knownColumns.push_back (std::make_shared<ColumnDataVector<std::uint16_t,BranchReader>> (branch->GetName()));
2278 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::int32_t>))
2279 {
2280 knownColumns.push_back (std::make_shared<ColumnDataVector<std::int32_t,BranchReader>> (branch->GetName()));
2281 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::uint32_t>))
2282 {
2283 knownColumns.push_back (std::make_shared<ColumnDataVector<std::uint32_t,BranchReader>> (branch->GetName()));
2284 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::int64_t>))
2285 {
2286 knownColumns.push_back (std::make_shared<ColumnDataVector<std::int64_t,BranchReader>> (branch->GetName()));
2287 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::uint64_t>))
2288 {
2289 knownColumns.push_back (std::make_shared<ColumnDataVector<std::uint64_t,BranchReader>> (branch->GetName()));
2290 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::vector<float>>))
2291 {
2292 knownColumns.push_back (std::make_shared<ColumnDataVectorVector<float,BranchReader>> (branch->GetName()));
2293 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::vector<std::int32_t>>))
2294 {
2295 knownColumns.push_back (std::make_shared<ColumnDataVectorVector<std::int32_t,BranchReader>> (branch->GetName()));
2296 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::vector<std::uint64_t>>))
2297 {
2298 knownColumns.push_back (std::make_shared<ColumnDataVectorVector<std::uint64_t,BranchReader>> (branch->GetName()));
2299 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::vector<std::vector<std::size_t>>>))
2300 {
2301 knownColumns.push_back (std::make_shared<ColumnDataVectorVectorVector<std::size_t,BranchReader>> (branch->GetName()));
2302 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::vector<std::vector<unsigned char>>>))
2303 {
2304 knownColumns.push_back (std::make_shared<ColumnDataVectorVectorVector<unsigned char,BranchReader>> (branch->GetName()));
2305 } else if (*branchClass->GetTypeInfo() == typeid(std::vector<std::string>))
2306 {
2307 knownColumns.push_back (std::make_shared<ColumnDataMetNames<BranchReader>> (branch->GetName()));
2308 }
2309 }
2310 }
2311 }
2312 }
2313 // This is a fallback for the case that we don't have an explicit
2314 // `samplingPattern` branch in our input file (i.e. an older file),
2315 // to allow us to still test tools needing it. This is likely not
2316 // something that actual users can do (they need the new files), but
2317 // for testing it seems like a reasonable workaround.
2318 knownColumns.push_back(std::make_shared<ColumnDataSamplingPattern<BranchReader>>("egammaClusters"));
2319
2320 // For branches that are element links they need to be explicitly
2321 // declared to have the correct xAOD type, correct split setting,
2322 // and correct linked containers.
2323
2324 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::CaloClusterContainer,BranchReader>>("AnalysisElectronsAuxDyn.caloClusterLinks"));
2325 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::TrackParticleContainer, BranchReader>>("AnalysisElectronsAuxDyn.trackParticleLinks"));
2326 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::CaloClusterContainer,BranchReader>>("AnalysisPhotonsAuxDyn.caloClusterLinks"));
2327 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::VertexContainer, BranchReader>>("AnalysisPhotonsAuxDyn.vertexLinks"));
2328 knownColumns.push_back(std::make_shared<ColumnDataVectorSplitLink<xAOD::TrackParticleContainer>>("AnalysisMuonsAuxDyn.inDetTrackParticleLink"));
2329 knownColumns.push_back(std::make_shared<ColumnDataVectorSplitLink<xAOD::TrackParticleContainer>>("AnalysisMuonsAuxDyn.combinedTrackParticleLink"));
2330 knownColumns.push_back(std::make_shared<ColumnDataVectorSplitLink<xAOD::TrackParticleContainer>>("AnalysisMuonsAuxDyn.extrapolatedMuonSpectrometerTrackParticleLink"));
2331 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::TrackParticleContainer, BranchReader>>("GSFConversionVerticesAuxDyn.trackParticleLinks"));
2332 knownColumns.push_back(std::make_shared<ColumnDataVectorSplitLink<xAOD::TrackParticleContainer>>("GSFTrackParticlesAuxDyn.originalTrackParticle"));
2333 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVariantLink<xAOD::IParticleContainer, BranchReader>>("AnalysisJetsAuxDyn.GhostTrack"));
2334 knownColumns.push_back(std::make_shared<ColumnDataVectorLink<xAOD::JetContainer, BranchReader>>("METAssoc_AnalysisMETAux.jetLink"));
2335 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVariantLink<xAOD::IParticleContainer, BranchReader>>("METAssoc_AnalysisMETAux.objectLinks"));
2336
2337 }else if (rntbackend)
2338 {
2339 std::unordered_map<std::string, ROOT::DescriptorId_t> fields;
2340 {
2341 const auto& desc = rntreader->GetDescriptor();
2342
2343 for (const auto& field : desc.GetTopLevelFields())
2344 {
2345 auto fieldName = field.GetFieldName();
2346 fields.emplace(desc.GetQualifiedFieldName(field.GetId()), field.GetId());
2347
2348 std::vector<ROOT::DescriptorId_t> subFieldIds{field.GetId()};
2349 while (!subFieldIds.empty())
2350 {
2351 const auto parentId = subFieldIds.back();
2352 auto parentname=desc.GetQualifiedFieldName(parentId);
2353 subFieldIds.pop_back();
2354
2355 for (const auto& subField : desc.GetFieldIterable(parentId))
2356 {
2357 auto subFieldName = desc.GetQualifiedFieldName(subField.GetId());
2358
2359 fields.emplace(desc.GetQualifiedFieldName(subField.GetId()), subField.GetId());
2360
2361 subFieldIds.push_back(subField.GetId());
2362 }
2363 }
2364 }
2365 }
2366
2367 const auto& desc = rntreader->GetDescriptor();
2368 for (const auto& [name, fieldId] : fields)
2369 {
2370 auto fieldName = desc.GetQualifiedFieldName(fieldId);
2371
2372 if (name.find("AuxDyn:") != std::string::npos ||
2373 name.find("Aux:") != std::string::npos)
2374 {
2375
2376 const auto& fieldDesc = desc.GetFieldDescriptor(fieldId);
2377 const std::string typeName = desc.GetTypeNameForComparison(fieldDesc);
2378 if (typeName == "std::int32_t" || typeName == "int")
2379 {
2380 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::int32_t, RNTFieldReader>>(name));
2381 } else if (typeName == "std::uint32_t" || typeName == "unsigned int")
2382 {
2383 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::uint32_t, RNTFieldReader>>(name));
2384 } else if (typeName == "std::uint64_t" || typeName == "unsigned long" || typeName == "unsigned long long")
2385 {
2386 knownColumns.push_back(std::make_shared<ColumnDataScalar<std::uint64_t, RNTFieldReader>>(name));
2387 } else if (typeName == "float")
2388 {
2389 knownColumns.push_back(std::make_shared<ColumnDataScalar<float, RNTFieldReader>>(name));
2390 } else if (typeName == "std::vector<float>")
2391 {
2392 knownColumns.push_back(std::make_shared<ColumnDataVector<float, RNTFieldReader>>(name));
2393 } else if (typeName == "std::vector<char>")
2394 {
2395 knownColumns.push_back(std::make_shared<ColumnDataVector<char, RNTFieldReader>>(name));
2396 } else if (typeName == "std::vector<std::int8_t>")
2397 {
2398 knownColumns.push_back(std::make_shared<ColumnDataVector<std::int8_t, RNTFieldReader>>(name));
2399 } else if (typeName == "std::vector<std::uint8_t>")
2400 {
2401 knownColumns.push_back(std::make_shared<ColumnDataVector<std::uint8_t, RNTFieldReader>>(name));
2402 } else if (typeName == "std::vector<std::int16_t>")
2403 {
2404 knownColumns.push_back(std::make_shared<ColumnDataVector<std::int16_t, RNTFieldReader>>(name));
2405 } else if (typeName == "std::vector<std::uint16_t>")
2406 {
2407 knownColumns.push_back(std::make_shared<ColumnDataVector<std::uint16_t, RNTFieldReader>>(name));
2408 } else if (typeName == "std::vector<std::int32_t>")
2409 {
2410 knownColumns.push_back(std::make_shared<ColumnDataVector<std::int32_t, RNTFieldReader>>(name));
2411 } else if (typeName == "std::vector<std::uint32_t>")
2412 {
2413 knownColumns.push_back(std::make_shared<ColumnDataVector<std::uint32_t, RNTFieldReader>>(name));
2414 } else if (typeName == "std::vector<std::int64_t>")
2415 {
2416 knownColumns.push_back(std::make_shared<ColumnDataVector<std::int64_t, RNTFieldReader>>(name));
2417 } else if (typeName == "std::vector<std::uint64_t>")
2418 {
2419 knownColumns.push_back(std::make_shared<ColumnDataVector<std::uint64_t, RNTFieldReader>>(name));
2420 } else if (typeName == "std::vector<std::vector<float>>")
2421 {
2422 knownColumns.push_back(std::make_shared<ColumnDataVectorVector<float, RNTFieldReader>>(name));
2423 } else if (typeName == "std::vector<std::vector<std::int32_t>>")
2424 {
2425 knownColumns.push_back(std::make_shared<ColumnDataVectorVector<std::int32_t, RNTFieldReader>>(name));
2426 } else if (typeName == "std::vector<std::vector<std::uint64_t>>")
2427 {
2428 knownColumns.push_back(std::make_shared<ColumnDataVectorVector<std::uint64_t, RNTFieldReader>>(name));
2429 } else if (typeName =="std::vector<std::vector<std::vector<std::size_t>>>")
2430 {
2431 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVector<std::size_t, RNTFieldReader>>(name));
2432 }else if (typeName =="std::vector<std::vector<std::vector<std::uint64_t>>>")
2433 {
2434 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVector<std::uint64_t, RNTFieldReader>>(name));
2435 }else if (typeName =="std::vector<std::vector<std::vector<std::uint8_t>>>")
2436 {
2437 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVector<std::uint8_t, RNTFieldReader>>(name));
2438 } else if (typeName =="std::vector<std::vector<std::vector<unsigned char>>>")
2439 {
2440 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVector<unsigned char, RNTFieldReader>>(name));
2441 } else if (typeName == "std::vector<std::string>")
2442 {
2443 knownColumns.push_back(std::make_shared<ColumnDataMetNames<RNTFieldReader>>(name));
2444 }
2445 }
2446 }
2447 knownColumns.push_back(std::make_shared<ColumnDataSamplingPattern<RNTFieldReader>>("egammaClusters"));
2448
2449 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::CaloClusterContainer,RNTFieldReader>>("AnalysisElectronsAuxDyn:caloClusterLinks"));
2450 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::TrackParticleContainer, RNTFieldReader>>("AnalysisElectronsAuxDyn:trackParticleLinks"));
2451 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::CaloClusterContainer,RNTFieldReader>>("AnalysisPhotonsAuxDyn:caloClusterLinks"));
2452 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::VertexContainer, RNTFieldReader>>("AnalysisPhotonsAuxDyn:vertexLinks"));
2453 knownColumns.push_back(std::make_shared<ColumnDataVectorRLink<xAOD::TrackParticleContainer,RNTFieldReader>>("AnalysisMuonsAuxDyn:inDetTrackParticleLink"));
2454 knownColumns.push_back(std::make_shared<ColumnDataVectorRLink<xAOD::TrackParticleContainer,RNTFieldReader>>("AnalysisMuonsAuxDyn:combinedTrackParticleLink"));
2455 knownColumns.push_back(std::make_shared<ColumnDataVectorRLink< xAOD::TrackParticleContainer, RNTFieldReader>>("AnalysisMuonsAuxDyn:extrapolatedMuonSpectrometerTrackParticleLink"));
2456 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorLink<xAOD::TrackParticleContainer, RNTFieldReader>>("GSFConversionVerticesAuxDyn:trackParticleLinks"));
2457 knownColumns.push_back(std::make_shared<ColumnDataVectorRLink<xAOD::TrackParticleContainer,RNTFieldReader>>("GSFTrackParticlesAuxDyn:originalTrackParticle"));
2458 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVariantLink<xAOD::IParticleContainer, RNTFieldReader>>("AnalysisJetsAuxDyn:GhostTrack"));
2459 knownColumns.push_back(std::make_shared<ColumnDataVectorLink<xAOD::JetContainer, RNTFieldReader>>("METAssoc_AnalysisMETAux:.jetLink"));
2460 knownColumns.push_back(std::make_shared<ColumnDataVectorVectorVariantLink<xAOD::IParticleContainer, RNTFieldReader>>("METAssoc_AnalysisMETAux:.objectLinks"));
2461
2462 }
2463
2464
2465 // For METMaker we need to preplace all of the MET terms that we
2466 // expect to be used, that's what this line does.
2467 std::vector<std::string> allMetTermNames;
2468 for (const auto& td : testDefinitions)
2469 {
2470 for (const auto& name : td.metTermNames)
2471 {
2472 if (std::find (allMetTermNames.begin(), allMetTermNames.end(), name) == allMetTermNames.end())
2473 allMetTermNames.push_back (name);
2474 }
2475 }
2476
2477
2478 if (tree)
2479 {
2480 if (!allMetTermNames.empty())
2481 knownColumns.push_back(std::make_shared<ColumnDataOutputMet<BranchReader>>("OutputMET",allMetTermNames));
2482
2483 // For METMaker we need various extra columns to run. This may need
2484 // some work to avoid, but would likey be worth it.
2485 knownColumns.push_back(std::make_shared<ColumnDataOutVector<std::uint16_t, BranchReader>>("AnalysisMuons.objectType", xAOD::Type::Muon));
2486 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, BranchReader>>("AnalysisMuons.m", ParticleConstants::muonMassInMeV));
2487 knownColumns.push_back(std::make_shared<ColumnDataOutVector<std::uint16_t, BranchReader>>("AnalysisJets.objectType", xAOD::Type::Jet));
2488
2489 // These are columns that represent variables that are normally held
2490 // by METAssociationHelper, or alternatively are decorated on the
2491 // MET terms (even though they are per object).
2492 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, BranchReader>>("AnalysisMuons.MetObjectWeight", 0));
2493 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, BranchReader>>("AnalysisJets.MetObjectWeight", 0));
2494 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, BranchReader>>("AnalysisJets.MetObjectWeightSoft", 0));
2495 knownColumns.push_back(std::make_shared<ColumnDataOutVector<MissingETBase::Types::bitmask_t,BranchReader>>("METAssoc_AnalysisMET.useObjectFlags", 0));
2496 } else if (rntbackend)
2497 {
2498 if (!allMetTermNames.empty())
2499 knownColumns.push_back(std::make_shared<ColumnDataOutputMet<BranchReader>>("OutputMET",allMetTermNames));
2500
2501 // For METMaker we need various extra columns to run. This may need
2502 // some work to avoid, but would likey be worth it.
2503 knownColumns.push_back(std::make_shared<ColumnDataOutVector<std::uint16_t, RNTFieldReader>>("AnalysisMuons.objectType", xAOD::Type::Muon));
2504 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, RNTFieldReader>>("AnalysisMuons.m", ParticleConstants::muonMassInMeV));
2505 knownColumns.push_back(std::make_shared<ColumnDataOutVector<std::uint16_t, RNTFieldReader>>("AnalysisJets.objectType", xAOD::Type::Jet));
2506
2507 // These are columns that represent variables that are normally held
2508 // by METAssociationHelper, or alternatively are decorated on the
2509 // MET terms (even though they are per object).
2510 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, RNTFieldReader>>("AnalysisMuons.MetObjectWeight", 0));
2511 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, RNTFieldReader>>("AnalysisJets.MetObjectWeight", 0));
2512 knownColumns.push_back(std::make_shared<ColumnDataOutVector<float, RNTFieldReader>>("AnalysisJets.MetObjectWeightSoft", 0));
2513 knownColumns.push_back(std::make_shared<ColumnDataOutVector<MissingETBase::Types::bitmask_t,RNTFieldReader>>("METAssoc_AnalysisMET.useObjectFlags", 0));
2514 }
2515 } // namespace columnar
constexpr double muonMassInMeV
the mass of the muon (in MeV)
@ Jet
The object is a jet.
Definition ObjectType.h:40
@ Muon
The object is a muon.
Definition ObjectType.h:48

Member Data Documentation

◆ file

std::unique_ptr<TFile> columnar::ColumnarPhysLiteTest::file

Definition at line 44 of file ColumnarPhysliteTest.h.

◆ inspector

std::unique_ptr<ROOT::Experimental::RNTupleInspector> columnar::ColumnarPhysLiteTest::inspector

Definition at line 47 of file ColumnarPhysliteTest.h.

◆ knownColumns

std::vector<std::shared_ptr<TestUtils::IColumnData> > columnar::ColumnarPhysLiteTest::knownColumns

Definition at line 49 of file ColumnarPhysliteTest.h.

◆ offsetColumns

std::unordered_map<std::string,const std::vector<ColumnarOffsetType>*> columnar::ColumnarPhysLiteTest::offsetColumns

Definition at line 51 of file ColumnarPhysliteTest.h.

◆ rntbackend

columnar::TestUtils::RNTupleBackend* columnar::ColumnarPhysLiteTest::rntbackend = nullptr

Definition at line 48 of file ColumnarPhysliteTest.h.

◆ rntreader

std::unique_ptr<ROOT::RNTupleReader> columnar::ColumnarPhysLiteTest::rntreader

Definition at line 46 of file ColumnarPhysliteTest.h.

◆ tree

TTree* columnar::ColumnarPhysLiteTest::tree = nullptr

Definition at line 45 of file ColumnarPhysliteTest.h.

◆ usedColumns

std::vector<std::shared_ptr<TestUtils::IColumnData> > columnar::ColumnarPhysLiteTest::usedColumns

Definition at line 50 of file ColumnarPhysliteTest.h.


The documentation for this struct was generated from the following files: