Jpp 20.0.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
JPinnacles.cc File Reference

Auxiliary program to filter acoustic fit results. More...

#include <iostream>
#include <iomanip>
#include <vector>
#include "TROOT.h"
#include "TFile.h"
#include "JAcoustics/JEvt.hh"
#include "JAcoustics/JSupport.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JFileRecorder.hh"
#include "JSupport/JMeta.hh"
#include "JTools/JRange.hh"
#include "Jeep/JContainer.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Auxiliary program to filter acoustic fit results.

Author
mdejong

Definition in file JPinnacles.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Auxiliary data structure for veto of acoustics events.

Apply veto to given acoustics event.

Parameters
evtacoustics event
Returns
true if vetoed; else false

Definition at line 27 of file JPinnacles.cc.

28{
29 using namespace std;
30 using namespace JPP;
31
32 /**
33 * Auxiliary data structure for veto of acoustics events.
34 */
35 struct veto_type :
36 public vector< JRange<double> >
37 {
38 /**
39 * Apply veto to given acoustics event.
40 *
41 * \param evt acoustics event
42 * \return true if vetoed; else false
43 */
44 bool operator()(const JEvt& evt) const
45 {
46 for (const_iterator i = this->begin(); i != this->end(); ++i) {
47 if ((*i)(evt.UNIXTimeStart) || (*i)(evt.UNIXTimeStop)) {
48 return true;
49 }
50 }
51
52 return false;
53 }
54 };
55
56
58 JLimit_t& numberOfEvents = inputFile.getLimit();
61 int debug;
62
63 try {
64
65 JParser<> zap("Auxiliary program to filter acoustic fit results.");
66
67 zap['f'] = make_field(inputFile, "input file (output of JKatoomba[.sh])");
68 zap['n'] = make_field(numberOfEvents) = JLimit::max();
69 zap['o'] = make_field(outputFile) = "pinnacles.root";
70 zap['V'] = make_field(veto, "veto acoustics events by UTC time ranges [s]");
71 zap['d'] = make_field(debug) = 2;
72
73 zap(argc, argv);
74 }
75 catch(const exception &error) {
76 FATAL(error.what() << endl);
77 }
78
79
80 outputFile.open();
81
82 outputFile.put(JMeta(argc, argv));
83
84 while (inputFile.hasNext()) {
85
86 STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
87
88 const JEvt* evt = inputFile.next();
89
90 if (!veto(*evt)) {
91 outputFile.put(*evt);
92 }
93 }
94
95 outputFile.close();
96}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
Object writing to file.
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Acoustic event fit.
double UNIXTimeStop
stop time
double UNIXTimeStart
start time
Auxiliary wrapper for I/O of container with optional comment (see JComment).
Definition JContainer.hh:42
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128
Auxiliary class for ROOT I/O of application specific meta data.
Definition JMeta.hh:72