Classes | Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | List of all members
AdcChannelDftPlotter Class Reference

#include <AdcChannelDftPlotter.h>

Inheritance diagram for AdcChannelDftPlotter:
AdcMultiChannelPlotter TpcDataTool AdcChannelTool

Classes

class  State
 
class  SubState
 

Public Types

using Index = unsigned int
 
using Name = std::string
 
- Public Types inherited from AdcMultiChannelPlotter
using Name = std::string
 
using NameVector = std::vector< Name >
 
using Index = unsigned int
 
using IndexVector = std::vector< Index >
 
using IndexSet = std::set< Index >
 
using AcdVector = std::vector< const AdcChannelData * >
 
using ChannelRangeMap = std::map< Name, IndexRange >
 
using ChannelGroupMap = std::map< Name, IndexRangeGroup >
 
using PadVector = std::vector< Pad >
 
- Public Types inherited from AdcChannelTool
using Index = unsigned int
 

Public Member Functions

 AdcChannelDftPlotter (fhicl::ParameterSet const &ps)
 
 ~AdcChannelDftPlotter () override
 
DataMap view (const AdcChannelData &acd) const override
 
int viewMapChannels (Name crn, const AcdVector &acds, TPadManipulator &man, Index ncr, Index icr) const override
 
int viewMapSummary (Index ilev, Name cgn, Name crn, TPadManipulator &man, Index ncr, Index icr) const override
 
bool updateWithView () const override
 
DataMap beginEvent (const DuneEventInfo &) const override
 
DataMap endEvent (const DuneEventInfo &) const override
 
- Public Member Functions inherited from AdcMultiChannelPlotter
 AdcMultiChannelPlotter (const fhicl::ParameterSet &ps, Name prefix="Plot")
 
 ~AdcMultiChannelPlotter () override
 
DataMap viewMap (const AdcChannelDataMap &acds) const override
 
void viewSummary (Index ilev) const
 
Index getLogLevel () const
 
Name getDataView () const
 
const NameVectorgetChannelRangeNames () const
 
const NameVectorgetChannelGroupNames () const
 
bool haveChannelRanges () const
 
bool haveChannelGroups () const
 
bool overlayGroups () const
 
Name getPlotName () const
 
Name getPlotSummaryName (Index ilev) const
 
Index getPlotSizeX () const
 
Index getPlotSizeY () const
 
Index getPlotSplitX () const
 
Index getPlotSplitY () const
 
const IndexRangeGroupgetChannelGroup (Name cgn) const
 
- Public Member Functions inherited from TpcDataTool
virtual DataMap updateTpcData (TpcData &) const
 
virtual DataMap viewTpcData (const TpcData &) const
 
virtual int forwardTpcData () const
 
- Public Member Functions inherited from AdcChannelTool
virtual ~AdcChannelTool ()=default
 
virtual DataMap update (AdcChannelData &) const
 
virtual DataMap updateMap (AdcChannelDataMap &acds) const
 
virtual bool viewWithUpdate () const
 
virtual DataMap close (const DataMap *dmin=nullptr)
 

Private Types

using IndexMap = std::map< Name, Index >
 
using HistMap = std::map< Name, TH1 * >
 

Private Member Functions

StategetState () const
 
SubStategetSubState (Index ilev) const
 
SubStategetJobState () const
 
SubStategetEventState () const
 
DataMap viewLocal (Name crn, const AcdVector &acds) const
 
int fillPad (DataMap &dm, TPadManipulator &man) const
 
bool skipChannel (const AdcChannelData &acd) const
 

Private Attributes

Name m_Variable
 
Index m_ChannelStatusFlag
 
Name m_ChannelSelection
 
float m_SampleFreq
 
double m_XMin
 
double m_XMax
 
float m_YMax
 
float m_YMinLog
 
Index m_NBinX
 
Name m_HistName
 
Name m_HistTitle
 
NameVector m_HistSummaryTitles
 
bool m_skipBad
 
bool m_skipNoisy
 
bool m_shiftFreq0
 
std::unique_ptr< TFormula > m_ptfsel
 
const AdcChannelStringToolm_adcStringBuilder
 
std::shared_ptr< Statem_pstate
 

Additional Inherited Members

- Static Public Member Functions inherited from AdcChannelTool
static int interfaceNotImplemented ()
 
- Protected Member Functions inherited from AdcMultiChannelPlotter
BaseStategetBaseState () const
 

Detailed Description

Definition at line 64 of file AdcChannelDftPlotter.h.

Member Typedef Documentation

using AdcChannelDftPlotter::HistMap = std::map<Name, TH1*>
private

Definition at line 117 of file AdcChannelDftPlotter.h.

using AdcChannelDftPlotter::Index = unsigned int

Definition at line 68 of file AdcChannelDftPlotter.h.

using AdcChannelDftPlotter::IndexMap = std::map<Name, Index>
private

Definition at line 116 of file AdcChannelDftPlotter.h.

Definition at line 69 of file AdcChannelDftPlotter.h.

Constructor & Destructor Documentation

AdcChannelDftPlotter::AdcChannelDftPlotter ( fhicl::ParameterSet const &  ps)

Definition at line 59 of file AdcChannelDftPlotter_tool.cc.

60 : AdcMultiChannelPlotter(ps, "Plot"),
61  m_Variable(ps.get<Name>("Variable")),
62  m_ChannelStatusFlag(ps.get<Index>("ChannelStatusFlag")),
63  m_ChannelSelection(ps.get<Name>("ChannelSelection")),
64  m_SampleFreq(ps.get<float>("SampleFreq")),
65  m_XMin(0.0),
66  m_XMax(0.0),
67  m_YMinLog(ps.get<float>("YMinLog")),
68  m_NBinX(0),
69  m_HistName(ps.get<Name>("HistName")),
70  m_HistTitle(ps.get<Name>("HistTitle")),
71  m_HistSummaryTitles(ps.get<NameVector>("HistSummaryTitles")),
72  m_pstate(new State) {
73  const string myname = "AdcChannelDftPlotter::ctor: ";
74  bool doMag = m_Variable == "magnitude";
75  bool doPha = m_Variable == "phase";
76  bool doPwr = m_Variable == "power";
77  bool doPwt = m_Variable == "power/tick";
78  // Check variable and get optional fields.
79  if ( doPwr || doPwt ) {
80  m_XMin = ps.get<double>("XMin");
81  m_XMax = ps.get<double>("XMax");
82  m_YMax = ps.get<float>("YMax");
83  m_NBinX = ps.get<Index>("NBinX");
84  } else if ( doMag ) {
85  m_YMax = ps.get<float>("YMax");
86  } else if ( doPha ) {
87  } else {
88  cout << myname << "Invalid Variable: " << m_Variable << endl;
89  throw art::Exception(art::errors::Configuration, "InvalidFclValue");
90  }
91  // Fetch renaming tools.
92  string snameBuilder = "adcStringBuilder";
95  if ( m_adcStringBuilder == nullptr ) {
96  cout << myname << "WARNING: AdcChannelStringTool not found: " << snameBuilder << endl;
97  }
98  // Derived config.
101  m_shiftFreq0 = (doPwr || doPwt) && (m_XMin >= m_XMax);
102  if ( m_ChannelSelection.size() ) {
103  m_ptfsel.reset(new TFormula("AdcChannelDftPlotter", m_ChannelSelection.c_str()));
104  }
105  // Display the configuration.
106  if ( getLogLevel() >= 1 ) {
107  cout << myname << "Configuration: " << endl;
108  cout << myname << " Variable: " << m_Variable << endl;
109  cout << myname << " ChannelStatusFlag: " << m_ChannelStatusFlag;
110  if ( m_skipBad ) {
111  if ( m_skipNoisy ) cout << " (skip bad and noisy)";
112  else cout << " (skip bad)";
113  } else if ( m_skipNoisy ) cout << " (skip noisy)";
114  cout << endl;
115  cout << myname << " ChannelSelection: " << m_ChannelSelection << endl;
116  cout << myname << " SampleFreq: " << m_SampleFreq << endl;
117  if ( doMag || doPwr || doPwt ) cout << myname << " NBinX: " << m_NBinX << endl;
118  if ( doPwr || doPwt ) {
119  cout << myname << " XMin: " << m_XMin << endl;
120  cout << myname << " XMax: " << m_XMax << endl;
121  cout << myname << " YMax: " << m_YMax << endl;
122  }
123  cout << myname << " YMinLog: " << m_YMinLog << endl;
124  cout << myname << " HistName: " << m_HistName << endl;
125  cout << myname << " HistTitle: " << m_HistTitle << endl;
126  cout << myname << " HistSummaryTitles: [";
127  bool first = true;
128  for ( Name name : m_HistSummaryTitles ) {
129  if ( first ) first = false;
130  else cout << ", ";
131  cout << name;
132  }
133  cout << "]" << endl;
134  }
135 }
static QCString name
Definition: declinfo.cpp:673
ChannelGroupService::Name Name
unsigned int Index
AdcMultiChannelPlotter(const fhicl::ParameterSet &ps, Name prefix="Plot")
std::unique_ptr< TFormula > m_ptfsel
static constexpr double ps
Definition: Units.h:99
std::shared_ptr< State > m_pstate
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::vector< string > NameVector
const AdcChannelStringTool * m_adcStringBuilder
static DuneToolManager * instance(std::string fclname="", int dbg=1)
T * getShared(std::string name)
QTextStream & endl(QTextStream &s)
AdcChannelDftPlotter::~AdcChannelDftPlotter ( )
override

Definition at line 139 of file AdcChannelDftPlotter_tool.cc.

139  {
140  const string myname = "AdcChannelDftPlotter::dtor: ";
141  if ( getLogLevel() >= 2 ) {
142  cout << myname << "Closing." << endl;
143  if ( getChannelRangeNames().size() ) {
144  cout << myname << "Channel ranges" << endl;
145  cout << myname << " CR name count nch/evt" << endl;
146  for ( Name crn : getChannelRangeNames() ) {
147  Index count = getJobState().count(crn);
148  double nchan = getJobState().nchan(crn);
149  cout << myname << setw(15) << crn << ":"
150  << setw(8) << count << setw(8) << nchan/count << endl;
151  }
152  } else {
153  cout << myname << "No channel ranges specified." << endl;
154  }
155  if ( getChannelRangeNames().size() ) {
156  for ( Name cgn : getChannelGroupNames() ) {
157  cout << myname << "Channel group " << cgn << endl;
158  cout << myname << " CR name count nch/evt" << endl;
159  const IndexRangeGroup& crg = getChannelGroup(cgn);
160  for ( const IndexRange& ran : crg.ranges ) {
161  Name crn = ran.name;
162  Index count = getJobState().count(crn);
163  double nchan = getJobState().nchan(crn);
164  cout << myname << setw(15) << crn << ":"
165  << setw(8) << count << setw(8) << nchan/count << endl;
166  }
167  }
168  } else {
169  cout << myname << "No channel groups specified." << endl;
170  }
171  }
172  viewSummary(0);
173 }
ChannelGroupService::Name Name
SubState & getJobState() const
unsigned int Index
const NameVector & getChannelRangeNames() const
Name name
Definition: IndexRange.h:32
void viewSummary(Index ilev) const
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
const IndexRangeGroup & getChannelGroup(Name cgn) const
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
RangeVector ranges
const NameVector & getChannelGroupNames() const
QTextStream & endl(QTextStream &s)

Member Function Documentation

DataMap AdcChannelDftPlotter::beginEvent ( const DuneEventInfo ) const
overridevirtual

Reimplemented from AdcChannelTool.

Definition at line 321 of file AdcChannelDftPlotter_tool.cc.

321  {
322  DataMap ret;
323  getEventState().clear();
324  return ret;
325 }
SubState & getEventState() const
DataMap AdcChannelDftPlotter::endEvent ( const DuneEventInfo ) const
overridevirtual

Reimplemented from AdcChannelTool.

Definition at line 329 of file AdcChannelDftPlotter_tool.cc.

329  {
330  DataMap ret;
331  viewSummary(1);
332  getEventState().clear();
333  return ret;
334 }
SubState & getEventState() const
void viewSummary(Index ilev) const
int AdcChannelDftPlotter::fillPad ( DataMap dm,
TPadManipulator man 
) const
private

Definition at line 529 of file AdcChannelDftPlotter_tool.cc.

529  {
530  const string myname = "AdcChannelDftPlotter::fillPad: ";
531  const bool dbg = false;
532  TGraph* pg = dm.getGraph("dftGraph");
533  TH1* ph = dm.getHist("dftHist");
534  float yValMax = dm.getFloat("dftYValMax");
535  bool doPha = m_Variable == "phase";
536  //bool doPwr = m_Variable == "power";
537  bool doPwt = m_Variable == "power/tick";
538  string dopt = dm.getString("dftDopt");
539  bool logy = false;
540  // Assign y limits.
541  double ymin = 0.0;
542  double ymax = 0.0;
543  if ( doPha ) {
544  ymax = 3.2;
545  ymin = -ymax;
546  } else {
547  ymin = 0.0;
548  if ( m_YMax > 0 ) ymax = m_YMax;
549  else if ( m_YMax < 0 && -m_YMax > yValMax ) ymax = -m_YMax;
550  else ymax = yValMax*1.02;
551  }
552  if ( m_YMinLog ) {
553  ymin = m_YMinLog;
554  logy = true;
555  }
556  double xmin = 0.0;
557  double xmax = 0.0;
558  Index ncr = dm.getInt("dftCRCount");
559  Index icr = 0;
560  bool manyCR = ncr > 1; // Does this plot have multiple CRs?
561  bool lastCR = true; // Is this the last CR on this plot?
562  if ( manyCR ) {
563  icr = dm.getInt("dftCRIndex");
564  lastCR = icr + 1 == ncr;
565  }
566  if ( pg != nullptr ) {
567  xmax = pg->GetXaxis()->GetXmax();
568  xmin = -0.02*xmax;
569  xmax *= 1.02;
570  if ( manyCR ) {
571  int icol = LineColors::color(icr, ncr);
572  pg->SetMarkerColor(icol);
573  }
574  man.add(pg, dopt);
575  } else if ( ph != nullptr ) {
576  if ( manyCR ) {
577  if ( icr > 0 ) dopt += " same";
578  int icol = LineColors::color(icr, ncr);
579  ph->SetLineColor(icol);
580  if ( dbg ) cout << myname << "DEBUG: Color[" << icr << "] = " << icol << ", dopt = " << dopt << endl;
581  }
582  man.add(ph, dopt);
583  } else {
584  cout << myname << "ERROR: Neither hist or graph is defined." << endl;
585  return 1;
586  }
587  if ( dbg ) cout << myname << "DEBUG: CR " << icr << "/" << ncr << endl;
588  // Build the descriptor strings.
589  string snevt;
590  if ( dm.haveInt("dftEventCount") ) {
591  ostringstream ssout;
592  ssout << "N_{ev} = " << dm.getInt("dftEventCount");
593  snevt = ssout.str();
594  }
595  string sncha; // # unique channels
596  string snven; // # view entries
597  {
598  ostringstream ssoutch;
599  ostringstream ssoutve;
600  ssoutch.precision(1);
601  ssoutve.precision(1);
602  if ( dm.haveFloat("dftChanPerEventCount") ) {
603  ssoutch << "N_{ch} = " << fixed << dm.getFloat("dftChanPerEventCount");
604  ssoutve << "N_{ve} = " << fixed << dm.getFloat("dftViewEntryPerEventCount");
605  } else {
607  const IntVector& allchans = dm.getIntVector("dftChannels");
608  Index nven = allchans.size();
609  Index ncha = 0;
610  for ( IntVector::const_iterator iven=allchans.begin(); iven!=allchans.end(); ++iven ) {
611  if ( find(allchans.begin(), iven, *iven) == iven ) ++ncha;
612  }
613  ssoutch << "N_{ch} = " << ncha;
614  ssoutve << "N_{ve} = " << nven;
615  // Display the view entry count if a view is defined
616  // or if the view entry and channl counts differ.
617  if ( getDataView().size() == 0 ) {
618  if ( nven != ncha ) {
619  cout << "ERROR: View entry count differs from channel count: "
620  << nven << " != " << ncha << "." << endl;
621  ssoutve << " !!!";
622  } else {
623  ssoutve.str("");
624  }
625  }
626  }
627  sncha = ssoutch.str();
628  snven = ssoutve.str();
629  }
630  string spow;
631  if ( doPwt ) {
632  ostringstream ssout;
633  double sum = ph->Integral(0, ph->GetNbinsX()+1);
634  ssout.precision(3);
635  ssout << "#sqrt{#Sigma} = " << fixed << setw(2) << sqrt(sum);
636  spow = ssout.str();
637  }
638  // If this is the last object added to the plot.
639  if ( lastCR ) {
640  if ( dbg ) cout << myname << "DEBUG: Closing plot." << endl;
641  man.addAxis();
642  if ( xmax > xmin ) man.setRangeX(xmin, xmax);
643  if ( ymax > ymin ) man.setRangeY(ymin, ymax);
644  if ( m_shiftFreq0 ) man.showUnderflow();
645  if ( logy ) man.setLogY();
646  if ( logy ) man.setGridY();
647  double textSize = 0.04;
648  int textFont = 42;
649  if ( ! manyCR ) {
650  double xlab = 0.70;
651  double ylab = 0.85;
652  double dylab = 1.2*textSize;
653  if ( spow.size() ) {
654  TLatex* ptxt = new TLatex(xlab, ylab, spow.c_str());
655  ptxt->SetNDC();
656  ptxt->SetTextFont(textFont);
657  ptxt->SetTextSize(textSize);
658  man.add(ptxt);
659  ylab -= dylab;
660  }
661  if ( sncha.size() ) {
662  TLatex* ptxt = new TLatex(xlab, ylab, sncha.c_str());
663  ptxt->SetNDC();
664  ptxt->SetTextFont(textFont);
665  ptxt->SetTextSize(textSize);
666  man.add(ptxt);
667  ylab -= dylab;
668  }
669  if ( snven.size() ) {
670  TLatex* ptxt = new TLatex(xlab, ylab, snven.c_str());
671  ptxt->SetNDC();
672  ptxt->SetTextFont(textFont);
673  ptxt->SetTextSize(textSize);
674  man.add(ptxt);
675  ylab -= dylab;
676  }
677  if ( snevt.size() ) {
678  TLatex* ptxt = new TLatex(xlab, ylab, snevt.c_str());
679  ptxt->SetNDC();
680  ptxt->SetTextFont(textFont);
681  ptxt->SetTextSize(textSize);
682  man.add(ptxt);
683  ylab -= dylab;
684  }
685  } else {
686  double xlab = 0.35;
687  double ylab = 0.85;
688  if ( snevt.size() ) {
689  TLatex* ptxt = new TLatex(xlab, ylab, snevt.c_str());
690  ptxt->SetNDC();
691  ptxt->SetTextFont(textFont);
692  ptxt->SetTextSize(textSize);
693  man.add(ptxt);
694  }
695  }
696  }
697  // Update legend.
698  if ( manyCR ) {
699  TObject* pobj = man.object();
700  Name lopt = pg == nullptr ? "l" : "p";
701  if ( icr == 0 ) {
702  double xlmin = 0.55;
703  double xlmax = 0.93;
704  double ylmax = 0.90;
705  double ylmin = ylmax - 0.05*(ncr+0.5);
706  if ( ylmin < 0.40 ) ylmin = 0.40;
707  man.addLegend(xlmin, ylmin, xlmax, ylmax);
708  } else {
709  pobj = man.object(icr);
710  }
711  TLegend* pleg = man.getLegend();
712  Name slab = dm.getString("dftCRLabel");
713  if ( pleg != nullptr ) {
714  pleg->SetMargin(0.1); // Fraction of box used for symbols
715  if ( spow.size() ) slab += " " + spow;
716  if ( sncha.size() ) slab += " " + sncha;
717  pleg->AddEntry(pobj, slab.c_str(), lopt.c_str());
718  }
719  }
720  return 0;
721 }
int setGridY(bool flag=true)
TObject * object() const
bool dbg
int add(unsigned int ipad, TObject *pobj, std::string sopt="", bool replace=false)
TLegend * addLegend(double x1, double y1, double x2, double y2)
ChannelGroupService::Name Name
bool haveInt(Name name) const
Definition: DataMap.h:207
const IntVector & getIntVector(Name name) const
Definition: DataMap.h:219
TGraph * getGraph(Name name) const
Definition: DataMap.h:225
intermediate_table::const_iterator const_iterator
static ColorType color(Index icolin, Index ncol=size())
Definition: LineColors.h:43
unsigned int Index
std::vector< int > IntVector
Definition: DataMap.h:50
int precision() const
Definition: qtextstream.h:259
String getString(Name name, String def="") const
Definition: DataMap.h:222
TLegend * getLegend() const
int setLogY(bool flag=true)
int showUnderflow(bool show=true)
int addAxis(bool flag=true)
Q_EXPORT QTSManip setw(int w)
Definition: qtextstream.h:331
Float getFloat(Name name, Float def=0.0) const
Definition: DataMap.h:220
std::vector< int > IntVector
Definition: fcldump.cxx:26
int getInt(Name name, int def=0) const
Definition: DataMap.h:218
bool haveFloat(Name name) const
Definition: DataMap.h:209
int setRangeX(double x1, double x2)
TH1 * getHist(Name name, TH1 *def=nullptr) const
Definition: DataMap.h:223
int setRangeY(double y1, double y2)
QTextStream & endl(QTextStream &s)
SubState& AdcChannelDftPlotter::getEventState ( ) const
inlineprivate

Definition at line 189 of file AdcChannelDftPlotter.h.

189 { return m_pstate->eventState(); };
std::shared_ptr< State > m_pstate
SubState& AdcChannelDftPlotter::getJobState ( ) const
inlineprivate

Definition at line 188 of file AdcChannelDftPlotter.h.

188 { return m_pstate->jobState(); };
std::shared_ptr< State > m_pstate
State& AdcChannelDftPlotter::getState ( ) const
inlineprivate

Definition at line 183 of file AdcChannelDftPlotter.h.

183 { return *m_pstate; };
std::shared_ptr< State > m_pstate
SubState& AdcChannelDftPlotter::getSubState ( Index  ilev) const
inlineprivate

Definition at line 184 of file AdcChannelDftPlotter.h.

184  {
185  if ( ilev > 1 ) abort();
186  return getState().m_ss[ilev];
187  }
bool AdcChannelDftPlotter::skipChannel ( const AdcChannelData acd) const
private

Definition at line 725 of file AdcChannelDftPlotter_tool.cc.

725  {
726  const string myname = "AdcChannelDftPlotter::skipChannel: ";
727  if ( ! m_ptfsel ) return false;
728  Index npar = m_ptfsel->GetNpar();
729  vector<double> pars(npar);
730  for ( Index ipar=0; ipar<npar; ++ipar ) {
731  string spar = m_ptfsel->GetParName(ipar);
732  if ( ! acd.hasAttribute(spar) ) {
733  cout << myname << "WARNING: Skipping run/event/channel " << acd.run() << "/" << acd.event()
734  << "/" << acd.channel() << " because of missing attribute " << spar << endl;
735  return true;
736  }
737  pars[ipar] = acd.getAttribute(spar);
738  }
739  double val = m_ptfsel->EvalPar(nullptr, &pars[0]);
740  return val == 0.0;
741 }
unsigned int Index
float getAttribute(Name mname, float def=0.0) const
std::unique_ptr< TFormula > m_ptfsel
AdcIndex run() const
AdcIndex event() const
bool hasAttribute(Name mname, float def=0.0) const
Channel channel() const
QTextStream & endl(QTextStream &s)
bool AdcChannelDftPlotter::updateWithView ( ) const
inlineoverridevirtual

Reimplemented from AdcChannelTool.

Definition at line 79 of file AdcChannelDftPlotter.h.

79 { return true; }
DataMap AdcChannelDftPlotter::view ( const AdcChannelData acd) const
overridevirtual

Reimplemented from AdcChannelTool.

Definition at line 304 of file AdcChannelDftPlotter_tool.cc.

304  {
305  const string myname = "AdcChannelDftPlotter::view: ";
306  AcdVector acds(1, &acd);
307  DataMap chret = viewLocal("", acds);
308  if ( getPlotName().size() ) {
310  TPadManipulator man;
311  fillPad(chret, man);
312  if ( getLogLevel() >= 3 ) cout << myname << "Printing " << pname << endl;
313  man.print(pname);
314  }
315  return chret;
316 }
AdcChannelDftPlotter::AcdVector AcdVector
DataMap viewLocal(Name crn, const AcdVector &acds) const
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
static std::string build(const AdcChannelStringTool *ptool, const AdcChannelData &acd, const DataMap &dm, std::string spat)
const AdcChannelStringTool * m_adcStringBuilder
int fillPad(DataMap &dm, TPadManipulator &man) const
int print(std::string fname, std::string spat="{,}")
QTextStream & endl(QTextStream &s)
DataMap AdcChannelDftPlotter::viewLocal ( Name  crn,
const AcdVector acds 
) const
private

Definition at line 338 of file AdcChannelDftPlotter_tool.cc.

338  {
339  const string myname = "AdcChannelDftPlotter::viewLocal: ";
340  if ( getLogLevel() >= 4 ) {
341  cout << myname << "Filling CR " << crn << " with ADC channel data size " << acds.size() << endl;
342  }
343  DataMap ret;
344  if ( acds.size() == 0 ) return ret;
345  const AdcChannelData* pacd = acds.front();
346  Index evt = pacd->event();
347  // Check if there have been any other views of this channel range for this event.
348  // For now, we implicity expect configurations that do not repeat ranges.
349  // Later we might cache results from an earlier attempt.
350  // For now, user can expect some plots to overcount events and maybe worse.
351  // Avoid this problem by ensuring configuration does not include any channel
352  // range more than once. This includes channel ranges in channel groups.
353  if ( getState().setEventChannelRange(evt, crn) ) {
354  cout << myname << "WARNING: Duplicate view of channel range " << crn
355  << " in event " << evt << endl;
356  }
357  if ( pacd == nullptr ) {
358  cout << myname << "ERROR: First channel has no data." << endl;
359  return ret;
360  }
361  const AdcChannelData& acd = *pacd;
362  bool doMag = m_Variable == "magnitude";
363  bool doPha = m_Variable == "phase";
364  bool doPwr = m_Variable == "power";
365  bool doPwt = m_Variable == "power/tick";
366  bool haveFreq = m_SampleFreq > 0.0;
367  if ( ! doMag && !doPha && !doPwr && !doPwt ) {
368  cout << myname << "ERROR: Invalid plot variable: " << m_Variable << endl;
369  return ret.setStatus(1);
370  }
371  Index nmag = acd.dftmags.size();
372  Index npha = acd.dftphases.size();
373  Index nsam = nmag + npha - 1;
374  if ( nmag == 0 ) {
375  cout << myname << "ERROR: DFT is not present." << endl;
376  return ret.setStatus(2);
377  }
378  if ( npha > nmag || nmag - npha > 1 ) {
379  cout << myname << "ERROR: DFT is not valid." << endl;
380  return ret.setStatus(3);
381  }
382  // Check consistency of input data.
383  Index nDataMissing = 0;
384  Index nBadMagCount = 0;
385  Index nBadPhaCount = 0;
386  for ( const AdcChannelData* pacde : acds ) {
387  if ( pacde == nullptr ) {
388  ++nDataMissing;
389  } else {
390  if ( pacde->dftmags.size() != nmag ) ++nBadMagCount;
391  if ( pacde->dftphases.size() != npha ) ++nBadPhaCount;
392  }
393  }
394  if ( nDataMissing ) cout << myname << "ERROR: Missing data channel count is " << nDataMissing << endl;
395  if ( nBadMagCount ) cout << myname << "ERROR: Inconsistent mag size channel count is " << nBadMagCount << endl;
396  if ( nBadPhaCount ) cout << myname << "ERROR: Inconsistent pha size channel count is " << nBadPhaCount << endl;
397  if ( nDataMissing || nBadMagCount || nBadPhaCount ) return ret.setStatus(4);
400  StringManipulator smanName(hname, false);
401  smanName.replace("%CRNAME%", crn);
402  smanName.replace("%VIEW%", getDataView());
403  StringManipulator smanTitl(htitl, false);
404  smanTitl.replace("%CRNAME%", crn);
405  smanTitl.replace("%VIEW%", getDataView());
406  //xx
407  //sman.replace("%CRNAME%", ran.name);
408  //sman.replace("%CRLABEL%", ran.label());
409  //sman.replace("%CRLABEL1%", ran.label(1));
410  //sman.replace("%CRLABEL2%", ran.label(2));
411  float pi = acos(-1.0);
412  double xFac = haveFreq ? m_SampleFreq/nsam : 1.0;
413  float yValMax = 0.0;
414  string dopt;
415  string xtitl = haveFreq ? "Frequency [kHz]" : "Frequency index";
416  if ( doMag || doPha ) {
417  if ( acds.size() != 1 ) {
418  cout << myname << "ERROR: " << (doMag ? "Magnitude" : "Phase")
419  << " may only be filled for a single channel." << endl;
420  return ret.setStatus(5);
421  }
422  string ytitl = "Phase";
423  if ( doMag ) {
424  ytitl = AdcChannelStringTool::build(m_adcStringBuilder, acd, "Amplitude% [SUNIT]%");
425  }
426  TGraph* pg = new TGraph;
427  pg->SetName(hname.c_str());
428  pg->SetTitle(htitl.c_str());
429  pg->SetMarkerStyle(2);
430  pg->SetMarkerColor(602);
431  //ph = new TH1F(hname.c_str(), htitl.c_str(), nbin, 0, nbin);
432  //ph->SetDirectory(nullptr);
433  //ph->SetLineWidth(2);
434  string xtitl = haveFreq ? "Frequency [kHz]" : "Frequency index";
435  pg->GetXaxis()->SetTitle(xtitl.c_str());
436  pg->GetYaxis()->SetTitle(ytitl.c_str());
437  for ( Index ipha=0; ipha<nmag; ++ipha ) {
438  float mag = acd.dftmags[ipha];
439  float pha = ipha < npha ? acd.dftphases[ipha] : 0.0;
440  if ( mag < 0 ) {
441  pha += ( pha > 0 ? -pi : pi);
442  mag = -mag;
443  }
444  float x = ipha*xFac;
445  float y = doPha ? pha : mag;
446  if ( y > yValMax ) yValMax = y;
447  pg->SetPoint(ipha, x, y);
448  }
449  ret.setGraph("dftGraph", pg);
450  ret.setString("dftDopt", "P");
451  } else if ( doPwr || doPwt ) {
452  // Build list of retained channels.
453  DataMap::IntVector dftChannels;
454  AcdVector keepAcds;
455  for ( const AdcChannelData* pacd : acds ) {
456  if ( m_ChannelStatusFlag ) {
457  if ( m_skipBad && pacd->channelStatus()==1 ) continue;
458  if ( m_skipNoisy && pacd->channelStatus()==2 ) continue;
459  }
460  if ( skipChannel(*pacd) ) continue;
461  dftChannels.push_back(pacd->channel());
462  keepAcds.push_back(pacd);
463  }
464  ret.setIntVector("dftChannels", dftChannels);
465  if ( getLogLevel() >= 4 ) {
466  Index ncha = channelCount(acds);
467  Index nchaKeep = channelCount(keepAcds);
468  cout << myname << " Retaining " << keepAcds.size() << " of " << acds.size() << " entries in "
469  << nchaKeep << " of " << ncha << " channels." << endl;
470  }
471  string ytitl = doPwr ? "Power" : "Power/tick";
472  if ( acd.sampleUnit.size() ) {
473  ytitl = AdcChannelStringTool::build(m_adcStringBuilder, acd, ytitl + " [(%SUNIT%)^{2}]");
474  }
475  double xmin = m_XMin;
476  double xmax = m_XMax;
477  Index nbin = m_NBinX;
478  // If min >= max, then show the full range.
479  // And for power plots shift so the zero frequency component is
480  // in the underflow bin and the highest is included in the last bin.
481  if ( xmin >= xmax ) {
482  xmin = 0.0;
483  xmax = (nmag-1)*xFac;
484  if ( m_shiftFreq0 ) {
485  double delx = 0.01*xFac;
486  xmin += delx;
487  xmax += delx;
488  }
489  }
490  if ( nbin == 0 ) {
491  nbin = (xmax - xmin)/xFac + 0.1;
492  }
493  TH1* ph = new TH1F(hname.c_str(), htitl.c_str(), nbin, xmin, xmax);
494  ph->SetDirectory(nullptr);
495  ph->SetLineWidth(2);
496  ph->GetXaxis()->SetTitle(xtitl.c_str());
497  ph->GetYaxis()->SetTitle(ytitl.c_str());
498  float pwrFac = 1.0/keepAcds.size();
499  if ( ! doPwr ) pwrFac /= nsam;
500  for ( Index imag=0; imag<nmag; ++imag ) {
501  float x = imag*xFac;
502  float y = 0.0;
503  for ( const AdcChannelData* pacde : keepAcds ) {
504  float mag = pacde->dftmags[imag];
505  if ( getLogLevel() >= 5 ) {
506  cout << myname << " Channel " << pacde->channel() << " sqrt(pwr[" << imag << "]): " << mag << endl;
507  }
508  y += pwrFac*mag*mag;
509  }
510  ph->Fill(x, y);
511  }
512  if ( ph->GetBinContent(nbin+1) && xmax > (nmag-1)*xFac ) {
513  cout << myname << "ERROR: Full range histogram has overflow." << endl;
514  }
515  for ( Index ibin=0; ibin<nbin; ++ibin ) {
516  double y = ph->GetBinContent(ibin);
517  if ( y > yValMax ) yValMax = y;
518  }
519  dopt = "hist";
520  ret.setHist("dftHist", ph, true);
521  ret.setString("dftDopt", "hist");
522  }
523  ret.setFloat("dftYValMax", yValMax);
524  return ret;
525 }
void setFloat(Name name, float val)
Definition: DataMap.h:133
DataMap & setStatus(int stat)
Definition: DataMap.h:130
AdcChannelDftPlotter::AcdVector AcdVector
bool skipChannel(const AdcChannelData &acd) const
void setHist(Name name, TH1 *ph, bool own=false)
Definition: DataMap.h:136
unsigned int Index
std::vector< int > IntVector
Definition: DataMap.h:50
void setIntVector(Name name, const IntVector &val)
Definition: DataMap.h:132
void setString(Name name, String val)
Definition: DataMap.h:135
static std::string build(const AdcChannelStringTool *ptool, const AdcChannelData &acd, const DataMap &dm, std::string spat)
AdcIndex event() const
AdcSignalVector dftphases
Index channelStatus() const
Channel channel() const
const AdcChannelStringTool * m_adcStringBuilder
list x
Definition: train.py:276
void setGraph(Name name, TGraph *pg)
Definition: DataMap.h:178
float pi
Definition: units.py:11
TCEvent evt
Definition: DataStructs.cxx:7
AdcSignalVector dftmags
QTextStream & endl(QTextStream &s)
int AdcChannelDftPlotter::viewMapChannels ( Name  crn,
const AcdVector acds,
TPadManipulator man,
Index  ncr,
Index  icr 
) const
overridevirtual

Implements AdcMultiChannelPlotter.

Definition at line 178 of file AdcChannelDftPlotter_tool.cc.

178  {
179  const string myname = "AdcChannelDftPlotter::viewMapChannels: ";
180  DataMap chret = viewLocal(crn, acds);
181  bool doState = true;
182  bool doFill = acds.size();
183  if ( doState ) {
184  ++getJobState().count(crn);
185  ++getEventState().count(crn);
186  Index jobCount = getJobState().count(crn);
187  Index evtCount = getEventState().count(crn);
188  bool doPwr = m_Variable == "power";
189  bool doPwt = m_Variable == "power/tick";
190  if ( doPwr || doPwt ) {
191  TH1* ph = chret.getHist("dftHist");
192  if ( ph != nullptr ) {
193  TH1*& phjob = getJobState().hist(crn);
194  if ( phjob == nullptr ) {
195  if ( jobCount == 1 ) {
196  phjob = dynamic_cast<TH1*>(ph->Clone());
197  phjob->SetDirectory(nullptr);
198  phjob->SetStats(0);
199  } else {
200  cout << myname << "ERROR: Hist missing for job count " << jobCount << endl;
201  }
202  } else {
203  if ( jobCount > 1 ) phjob->Add(ph);
204  else cout << myname << "ERROR: Hist existing for job count " << jobCount << endl;
205  }
206  TH1*& phevt = getEventState().hist(crn);
207  if ( phevt == nullptr ) {
208  if ( evtCount == 1 ) {
209  phevt = dynamic_cast<TH1*>(ph->Clone());
210  phevt->SetDirectory(nullptr);
211  phevt->SetStats(0);
212  } else {
213  cout << myname << "ERROR: Hist missing for event count " << evtCount << endl;
214  }
215  } else {
216  if ( evtCount > 1 ) phevt->Add(ph);
217  else cout << myname << "ERROR: Hist existing for event count " << evtCount << endl;
218  }
219  const IntVector& allchans = chret.getIntVector("dftChannels");
220  Index nven = allchans.size();
221  Index ncha = channelCount(allchans);
222  getJobState().nchan(crn) += ncha;
223  getJobState().nviewentry(crn) += nven;
224  getEventState().nchan(crn) += ncha;
225  getEventState().nviewentry(crn) += nven;
226  } else {
227  doFill = false; // No data in CR
228  }
229  }
230  }
231  chret.setString("dftCRLabel", crn);
232  chret.setInt("dftCRCount", 1);
233  chret.setInt("dftCRIndex", icr);
234  if ( doFill ) fillPad(chret, man);
235  return 0;
236 }
SubState & getJobState() const
const IntVector & getIntVector(Name name) const
Definition: DataMap.h:219
SubState & getEventState() const
DataMap viewLocal(Name crn, const AcdVector &acds) const
unsigned int Index
void setString(Name name, String val)
Definition: DataMap.h:135
void setInt(Name name, int val)
Definition: DataMap.h:131
std::vector< int > IntVector
Definition: fcldump.cxx:26
TH1 * getHist(Name name, TH1 *def=nullptr) const
Definition: DataMap.h:223
int fillPad(DataMap &dm, TPadManipulator &man) const
QTextStream & endl(QTextStream &s)
int AdcChannelDftPlotter::viewMapSummary ( Index  ilev,
Name  cgn,
Name  crn,
TPadManipulator man,
Index  ncr,
Index  icr 
) const
overridevirtual

Implements AdcMultiChannelPlotter.

Definition at line 241 of file AdcChannelDftPlotter_tool.cc.

241  {
242  const string myname = "AdcChannelDftPlotter::viewMapSummary: ";
243  if ( getLogLevel() >= 2 ) {
244  cout << myname << "Processing " << cgn << "/" << crn << " (" << icr << "/" << ncr << ")" << endl;
245  }
246  if ( ilev > 1 ) {
247  cout << myname << "ERROR: Invalid level: " << ilev << endl;
248  return 12;
249  }
250  if ( icr >= ncr ) {
251  cout << myname << "ERROR: Too many plots: " << icr << " >= " << ncr << endl;
252  return 11;
253  }
254  SubState& levState = getSubState(ilev);
255  Index count = levState.count(crn);
256  Index nchanTot = levState.nchan(crn);
257  float nchanEvt = count > 0 ? double(nchanTot)/count : 0.0;
258  Index nvenTot = levState.nviewentry(crn);
259  float nvenEvt = count > 0 ? double(nvenTot)/count : 0.0;
260  TH1* ph = nullptr;
261  TH1* phin = levState.hist(crn);
262  bool doFill = false;
263  if ( phin != nullptr ) {
264  if ( count == 0 ) return 12;
265  ph = (phin == nullptr) ? nullptr : dynamic_cast<TH1*>(phin->Clone());
266  if ( ph == nullptr ) return 13;
267  ph->SetDirectory(nullptr);
268  if ( m_HistSummaryTitles.size() ) {
269  Name htitl = ilev < m_HistSummaryTitles.size()
270  ? m_HistSummaryTitles[ilev]
271  : m_HistSummaryTitles.back();
272  StringManipulator smanTitl(htitl, false);
273  Name cglab = getChannelGroup(cgn).label();
274  if ( cglab.size() == 0 ) cglab = cgn;
275  smanTitl.replace("%CGNAME%", cgn);
276  smanTitl.replace("%CGLABEL%", cglab);
277  smanTitl.replace("%CRNAME%", crn);
278  smanTitl.replace("%RUN%", getBaseState().run());
279  smanTitl.replace("%EVENT%", getBaseState().event);
280  smanTitl.replace("%VIEW%", getDataView());
281  ph->SetTitle(htitl.c_str());
282  }
283  double fac = 1.0/count;
284  ph->Scale(fac);
285  doFill = true;
286  }
287  DataMap dm;
288  dm.setHist("dftHist", ph, true);
289  dm.setInt("dftEventCount", count);
290  dm.setFloat("dftChanPerEventCount", nchanEvt);
291  dm.setFloat("dftViewEntryPerEventCount", nvenEvt);
292  dm.setString("dftDopt", "hist");
293  dm.setString("dftCRLabel", crn);
294  dm.setInt("dftCRCount", ncr);
295  dm.setInt("dftCRIndex", icr);
296  if ( doFill ) fillPad(dm, man);
297  //man.add(ph, "hist");
298  //delete ph;
299  return 0;
300 }
SubState & getSubState(Index ilev) const
void setFloat(Name name, float val)
Definition: DataMap.h:133
BaseState & getBaseState() const
Name label(Index ilab=0) const
ChannelGroupService::Name Name
void setHist(Name name, TH1 *ph, bool own=false)
Definition: DataMap.h:136
unsigned int Index
const IndexRangeGroup & getChannelGroup(Name cgn) const
void setString(Name name, String val)
Definition: DataMap.h:135
void setInt(Name name, int val)
Definition: DataMap.h:131
int fillPad(DataMap &dm, TPadManipulator &man) const
QTextStream & endl(QTextStream &s)
Event finding and building.

Member Data Documentation

const AdcChannelStringTool* AdcChannelDftPlotter::m_adcStringBuilder
private

Definition at line 106 of file AdcChannelDftPlotter.h.

Name AdcChannelDftPlotter::m_ChannelSelection
private

Definition at line 88 of file AdcChannelDftPlotter.h.

Index AdcChannelDftPlotter::m_ChannelStatusFlag
private

Definition at line 87 of file AdcChannelDftPlotter.h.

Name AdcChannelDftPlotter::m_HistName
private

Definition at line 95 of file AdcChannelDftPlotter.h.

NameVector AdcChannelDftPlotter::m_HistSummaryTitles
private

Definition at line 97 of file AdcChannelDftPlotter.h.

Name AdcChannelDftPlotter::m_HistTitle
private

Definition at line 96 of file AdcChannelDftPlotter.h.

Index AdcChannelDftPlotter::m_NBinX
private

Definition at line 94 of file AdcChannelDftPlotter.h.

std::shared_ptr<State> AdcChannelDftPlotter::m_pstate
private

Definition at line 182 of file AdcChannelDftPlotter.h.

std::unique_ptr<TFormula> AdcChannelDftPlotter::m_ptfsel
private

Definition at line 103 of file AdcChannelDftPlotter.h.

float AdcChannelDftPlotter::m_SampleFreq
private

Definition at line 89 of file AdcChannelDftPlotter.h.

bool AdcChannelDftPlotter::m_shiftFreq0
private

Definition at line 102 of file AdcChannelDftPlotter.h.

bool AdcChannelDftPlotter::m_skipBad
private

Definition at line 100 of file AdcChannelDftPlotter.h.

bool AdcChannelDftPlotter::m_skipNoisy
private

Definition at line 101 of file AdcChannelDftPlotter.h.

Name AdcChannelDftPlotter::m_Variable
private

Definition at line 86 of file AdcChannelDftPlotter.h.

double AdcChannelDftPlotter::m_XMax
private

Definition at line 91 of file AdcChannelDftPlotter.h.

double AdcChannelDftPlotter::m_XMin
private

Definition at line 90 of file AdcChannelDftPlotter.h.

float AdcChannelDftPlotter::m_YMax
private

Definition at line 92 of file AdcChannelDftPlotter.h.

float AdcChannelDftPlotter::m_YMinLog
private

Definition at line 93 of file AdcChannelDftPlotter.h.


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