Skip to content

Commit

Permalink
Merge pull request #831 from JeffersonLab/aaust_FCAL_hits
Browse files Browse the repository at this point in the history
Aaust fcal hits
  • Loading branch information
sdobbs authored Aug 22, 2024
2 parents 6eb0575 + 66725c6 commit 8e9304a
Show file tree
Hide file tree
Showing 5 changed files with 77 additions and 0 deletions.
47 changes: 47 additions & 0 deletions src/libraries/HDDM/DEventSourceREST.cc
Original file line number Diff line number Diff line change
Expand Up @@ -398,6 +398,10 @@ jerror_t DEventSourceREST::GetObjects(JEvent &event, JFactory_base *factory)
return Extract_DFMWPCHit(record,
dynamic_cast<JFactory<DFMWPCHit>*>(factory), locEventLoop);
}
if (dataClassName =="DFCALHit") {
return Extract_DFCALHit(record,
dynamic_cast<JFactory<DFCALHit>*>(factory), locEventLoop);
}
if (dataClassName =="DDetectorMatches") {
return Extract_DDetectorMatches(locEventLoop, record,
dynamic_cast<JFactory<DDetectorMatches>*>(factory));
Expand Down Expand Up @@ -2029,6 +2033,49 @@ jerror_t DEventSourceREST::Extract_DFMWPCHit(hddm_r::HDDM *record,
return NOERROR;
}

//-----------------------
// Extract_DFCALHit
//-----------------------
jerror_t DEventSourceREST::Extract_DFCALHit(hddm_r::HDDM *record,
JFactory<DFCALHit>* factory, JEventLoop* locEventLoop)
{
/// Copies the data from the fcal hit hddm record. This is
/// call from JEventSourceREST::GetObjects. If factory is NULL, this
/// returns OBJECT_NOT_AVAILABLE immediately.

if (factory==NULL) {
return OBJECT_NOT_AVAILABLE;
}
string tag = (factory->Tag())? factory->Tag() : "";

vector<DFCALHit*> data;

// loop over fmwpc hit records
const hddm_r::FcalHitList &hits =
record->getFcalHits();
hddm_r::FcalHitList::iterator iter;
for (iter = hits.begin(); iter != hits.end(); ++iter) {
if (iter->getJtag() != tag)
continue;

DFCALHit *hit = new DFCALHit();
hit->row = iter->getRow();
hit->column = iter->getColumn();
hit->x = iter->getX();
hit->y = iter->getY();
hit->E = iter->getE();
hit->t = iter->getT();
hit->intOverPeak = iter->getIntOverPeak();

data.push_back(hit);
}

// Copy into factory
factory->CopyTo(data);

return NOERROR;
}

//----------------------------
// Extract_DEventHitStatistics
//----------------------------
Expand Down
3 changes: 3 additions & 0 deletions src/libraries/HDDM/DEventSourceREST.h
Original file line number Diff line number Diff line change
Expand Up @@ -28,6 +28,7 @@
#include <TRACKING/DTrackTimeBased.h>
#include <FCAL/DFCALShower.h>
#include <FCAL/DFCALShower_factory.h>
#include <FCAL/DFCALHit.h>
#include <CCAL/DCCALShower.h>
#include <BCAL/DBCALShower.h>
#include <BCAL/DBCALShower_factory_IU.h>
Expand Down Expand Up @@ -101,6 +102,8 @@ class DEventSourceREST:public JEventSource
JFactory<DDIRCPmtHit>* factory, JEventLoop* locEventLoop);
jerror_t Extract_DFMWPCHit(hddm_r::HDDM *record,
JFactory<DFMWPCHit>* factory, JEventLoop* locEventLoop);
jerror_t Extract_DFCALHit(hddm_r::HDDM *record,
JFactory<DFCALHit>* factory, JEventLoop* locEventLoop);
jerror_t Extract_DEventHitStatistics(hddm_r::HDDM *record,
JFactory<DEventHitStatistics> *factory);

Expand Down
23 changes: 23 additions & 0 deletions src/libraries/HDDM/DEventWriterREST.cc
Original file line number Diff line number Diff line change
Expand Up @@ -57,6 +57,9 @@ DEventWriterREST::DEventWriterREST(JEventLoop* locEventLoop, string locOutputFil
ADD_FCAL_DATA_FOR_CPP=false;
gPARMS->SetDefaultParameter("PID:ADD_FCAL_DATA_FOR_CPP",ADD_FCAL_DATA_FOR_CPP);

REST_WRITE_FCAL_HITS = false;
gPARMS->SetDefaultParameter("REST:WRITE_FCAL_HITS", REST_WRITE_FCAL_HITS);


CCDB_CONTEXT_STRING = "";
// if we can get the calibration context from the DANA interface, then save this as well
Expand Down Expand Up @@ -88,6 +91,11 @@ bool DEventWriterREST::Write_RESTEvent(JEventLoop* locEventLoop, string locOutpu
std::vector<const DFCALShower*> fcalshowers;
locEventLoop->Get(fcalshowers);

std::vector<const DFCALHit*> fcalhits;
if(REST_WRITE_FCAL_HITS) {
locEventLoop->Get(fcalhits);
}

std::vector<const DBCALShower*> bcalshowers;
locEventLoop->Get(bcalshowers);

Expand Down Expand Up @@ -428,6 +436,21 @@ bool DEventWriterREST::Write_RESTEvent(JEventLoop* locEventLoop, string locOutpu
}
}

if(REST_WRITE_FCAL_HITS) {
// push any DFCALHit objects to the output record
for (size_t i=0; i < fcalhits.size(); i++)
{
hddm_r::FcalHitList hit = res().addFcalHits(1);
hit().setRow(fcalhits[i]->row);
hit().setColumn(fcalhits[i]->column);
hit().setX(fcalhits[i]->x);
hit().setY(fcalhits[i]->y);
hit().setE(fcalhits[i]->E);
hit().setT(fcalhits[i]->t);
hit().setIntOverPeak(fcalhits[i]->intOverPeak);
}
}

// push any DTrackTimeBased objects to the output record
for (size_t i=0; i < tracks.size(); ++i)
{
Expand Down
2 changes: 2 additions & 0 deletions src/libraries/HDDM/DEventWriterREST.h
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,7 @@
#include "PID/DBeamPhoton.h"
#include "TRACKING/DMCThrown.h"
#include "FCAL/DFCALShower.h"
#include "FCAL/DFCALHit.h"
#include "CCAL/DCCALShower.h"
#include "PID/DNeutralShower.h"
#include <PID/DDetectorMatches.h>
Expand Down Expand Up @@ -62,6 +63,7 @@ class DEventWriterREST : public JObject
bool REST_WRITE_CCAL_SHOWERS;
bool REST_WRITE_TRACK_EXIT_PARAMS;
bool ADD_FCAL_DATA_FOR_CPP;
bool REST_WRITE_FCAL_HITS;

// metadata to save in the REST file
// these should be consistent during program execution
Expand Down
2 changes: 2 additions & 0 deletions src/libraries/HDDM/rest.xml
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,8 @@
E1E9="float" E9E25="float"/>
<fcalShowerNBlocks minOccurs="0" numBlocks="int"/>
</fcalShower>
<fcalHit maxOccurs="unbounded" minOccurs="0" jtag="string" row="int" column="int"
x="float" y="float" E="float" Eunit="GeV" t="float" tunit="ns" intOverPeak="float"/>
<bcalShower minOccurs="0" maxOccurs="unbounded" jtag="string"
x="float" y="float" z="float" t="float" E="float"
xerr="float" yerr="float" zerr="float" terr="float" Eerr="float"
Expand Down

0 comments on commit 8e9304a

Please sign in to comment.