1 #ifndef __JSUPPORT__JTREESCANNER__
2 #define __JSUPPORT__JTREESCANNER__
27 namespace JPP {
using namespace JSUPPORT; }
51 template<
class JClass_t>
53 public JPointer< JChainReader<JClass_t> >
76 gErrorIgnoreLevel = kError;
90 template<
class JClass_t = JNullType,
class JEvaluator_t = JNullType>
112 template<
class JDerived_t,
class JBase_t>
174 if (counter < this->getLowerLimit()) {
175 skip(this->getLowerLimit() - counter);
178 return (counter < this->getEntries() &&
179 counter < this->getUpperLimit());
190 if (this->hasNext()) {
192 this->
get()->GetEvent(counter++);
194 ps.reset(this->
get()->getAddress());
225 using namespace JROOT;
227 this->setLimit(limit);
231 this->
get()->Reset();
235 TFile* file = TFile::Open(i->c_str());
239 if (file->GetListOfKeys()->Contains(this->
get()->getTreeName())) {
240 this->
get()->Add(i->c_str());
249 actionAtFileOpen<JDerived_t>(this->
get()->GetCurrentFile());
250 actionAtFileOpen<JBase_t> (this->
get()->GetCurrentFile());
261 return this->
get()->GetEntries();
273 if (index >= 0 && index < this->getEntries()) {
275 this->
get()->GetEvent(index);
277 return this->
get()->getAddress();
302 return JDerived_t::Class_Name();
317 template<
class JClass_t>
319 public JTreeScanner<JAssertConversion<JClass_t, JClass_t>, JNullType>
362 template<
class JDerived_t,
class JBase_t,
class JEvaluator_t>
365 public JTreeScanner<JAssertConversion<JDerived_t, JBase_t>, JNullType>
387 const double __value) :
402 return first.value < second.value;
413 friend inline bool operator<(
const JEntry_t& entry,
const double value)
415 return entry.value < value;
456 const JEvaluator_t& evaluator = JEvaluator_t()) :
472 const JEvaluator_t& evaluator = JEvaluator_t()) :
486 if (this->hasNext()) {
488 this->
get()->GetEvent(
queue[this->counter++].index);
490 ps.reset(this->
get()->getAddress());
513 setBranchStatus(this->
get()->GetBranch(this->
get()->getBranchName()),
BRANCH_PATTERN,
false);
515 queue.resize(this->getEntries());
517 typename queue_type::iterator out =
queue.begin();
519 for (Long64_t i = 0, n0 = 0; i != this->getEntries(); ++i, ++out) {
521 const Long64_t n1 = (100 * (i + 1)) / this->getEntries();
525 STATUS(left << setw(24) << this->
get()->GetName() << right <<
' ' << setw(3) << n1 <<
"%\r");
DEBUG(endl);
530 this->
get()->GetEvent(i);
532 const data_type* p = this->
get()->getAddress();
538 this->
get()->SetBranchStatus(
"*", 1);
552 if (index >= 0 && index < (Long64_t)
queue.size()) {
554 this->
get()->GetEvent(
queue[index].index);
556 return this->
get()->getAddress();
571 virtual Long64_t
find(
const double value)
const override
575 if (!
queue.empty()) {
577 typename queue_type::const_iterator p = lower_bound(
queue.begin(),
queue.end(), value);
579 if (p ==
queue.end()) {
581 return queue.size() - 1;
583 }
else if (p ==
queue.begin()) {
589 typename queue_type::const_iterator q = p--;
591 if (value - p->value < q->value - value)
615 if (branch != NULL) {
617 TObjArray* array = branch->GetListOfBranches();
619 for (Int_t i = 0; i != array->GetEntries(); ++i) {
621 TBranch* p = (TBranch*) array->At(i);
623 if (p != NULL &&
string(p->GetName()).find(pattern) != string::npos) {
625 if (p->GetSplitLevel() == 0) {
627 NOTICE(
"Set status of branch " << p->GetName() <<
" to " << status << endl);
629 p->SetStatus(status);
633 setBranchStatus(p, pattern, status);
650 template<
class JClass_t,
class JEvaluator_t>
652 public JTreeScanner<JAssertConversion<JClass_t, JClass_t>, JEvaluator_t>
668 const JEvaluator_t& evaluator = JEvaluator_t()) :
682 const JEvaluator_t& evaluator = JEvaluator_t()) :
virtual const pointer_type & next() override
Get next element.
virtual void configure(const JMultipleFileScanner_t &file_list, const JLimit &limit) override
Configure.
virtual Long64_t find(const double value) const override
Find index of element that is closest in value to given value.
JTreeScanner(const JMultipleFileScanner_t &file_list, const JLimit &limit=JLimit())
Constructor.
JTreeScanner(const JMultipleFileScanner_t &file_list, const JLimit &limit, const JEvaluator_t &evaluator=JEvaluator_t())
Constructor.
JEntry_t(const Long64_t __index, const double __value)
Constructor.
virtual bool hasNext() override
Check availability of next element.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
unsigned int skip_type
Type definition for number of objects to skip.
JTreeScanner(const JMultipleFileScanner_t &file_list, const JEvaluator_t &evaluator=JEvaluator_t())
Constructor.
then echo Submitting reweighting and histogram comparison jobs to nikhef stbc batch queue
static void setBranchStatus(TBranch *branch, const char *pattern, const bool status)
Set status of branch.
JTreeScannerInterface< JBase_t >::pointer_type pointer_type
JTreeScanner(const JMultipleFileScanner_t &file_list, const JLimit &limit, const JEvaluator_t &evaluator=JEvaluator_t())
Constructor.
Type definition for counter for ROOT TTree and auxiliary methods.
JTreeScanner(const JMultipleFileScanner_t &file_list, const JEvaluator_t &evaluator=JEvaluator_t())
Constructor.
Auxialiary class to assert type conversion.
Long64_t counter_type
Type definition for counter.
friend bool operator<(const JEntry_t &entry, const double value)
Comparison between TChain entry and value.
virtual void rewind() override
Rewind.
JTreeScanner()
Default constructor.
bool is_valid() const
Check validity of pointer.
Template definition for direct access of elements in ROOT TChain.
JTreeScanner()
Default constructor.
virtual skip_type skip(const skip_type ns) override
Skip items.
Auxiliary interface for direct access of elements in ROOT TChain.
then echo The file $DIR KM3NeT_00000001_00000000 root already please rename or remove it first
Auxiliary class for template TChain reading.
virtual data_type * getEntry(Long64_t index) override
Get entry at given index.
virtual void configure(const JMultipleFileScanner_t &file_list, const JLimit &limit) override
Configure.
virtual Long64_t getEntries() const override
Get number of entries.
Scanning of objects from a single file according a format that follows from the extension of each fil...
virtual data_type * getEntry(Long64_t index) override
Get entry at given index.
friend bool operator<(const JEntry_t &first, const JEntry_t &second)
Comparison between two TChain entries.
std::ostream & rewind(std::ostream &out)
Rewind character.
virtual const pointer_type & next() override
Get next element.
Auxiliary class for defining the range of iterations of objects.
JTreeScannerInterface< JBase_t, JEvaluator_t >::pointer_type pointer_type
TChain reading for template data type.
JTreeScanner()
Default constructor.
Long64_t index
index in TChain
std::vector< JEntry_t > queue_type
Type definition of internal queue for ordering the elements in the TChain.
Auxiliary class for no type definition.
Template implementation of class that holds pointer to object(s).
Data type dependent action methods for customised ROOT version management.
virtual counter_type getCounter() const override
Get internal counter.
virtual const char * getClassName() const override
Get actual class name.
General purpose messaging.
counter_type advance(counter_type &counter, const counter_type value, const counter_type limit=std::numeric_limits< counter_type >::max())
Advance counter.
JTreeScanner(const JTreeScanner &input)
Copy constructor.
static const char *const BRANCH_PATTERN
Pattern match for names of sub-branches that will not be read when ordering elements in a TTree...
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Base class for JTreeScanner.
JEntry_t()
Default constructor.
Auxiliary base class for list of file names.
JTreeScanner(const JTreeScanner &input)
Copy constructor.
General purpose class for object reading from a list of file names.
JTreeScanner(const JMultipleFileScanner_t &file_list, const JLimit &limit=JLimit())
Constructor.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
double value
corresponding value
JTreeScanner()
Default constructor.
JTreeScanner_t()
Default constructor.
Auxiliary class for handling debug parameter within a class.