23 #include <Geant4/G4SystemOfUnits.hh>
24 #include <Geant4/Randomize.hh>
27 #include <TSQLServer.h>
28 #include <TSQLResult.h>
32 #include <boost/lexical_cast.hpp>
38 if(
Verbosity() > 2) std::cout <<
"SQDigitizer: Init " << detIDByName.size() << std::endl;
43 if(!enableDC1 && (i >= 7 && i <= 12))
continue;
47 detIDByName[detName] = i;
59 digitize_secondaries(false)
70 if(
Verbosity() > 2) std::cout <<
"SQDigitizer: InitRun " << detIDByName.size() << std::endl;
79 std::cerr <<
Name() <<
" DST Node missing, abort." << std::endl;
84 hitContainerByName.clear();
85 for(
auto it = detIDByName.begin(); it != detIDByName.end(); ++it)
88 std::string detName = it->first;
89 std::string g4hitNodeName =
"G4HIT_" + detName;
92 std::cout <<
Name() <<
": booking input G4HIT node " << g4hitNodeName << std::endl;
95 PHG4HitContainer* hits = findNode::getClass<PHG4HitContainer>(topNode, g4hitNodeName.c_str());
98 if(
Verbosity() > 2) std::cout <<
Name() <<
": Could not locate g4 hit node " << g4hitNodeName << std::endl;
103 hitContainerByName[detName] = hits;
108 if(hitContainerByName.empty())
110 std::cerr <<
Name() <<
": no g4hit node foud, abort " << std::endl;
115 digits = findNode::getClass<SQHitVector>(topNode,
"SQHitVector");
118 if(
Verbosity() > 2) std::cout <<
Name() <<
": booking output node SQHitVector" << std::endl;
130 for(
auto it = hitContainerByName.begin(); it != hitContainerByName.end(); ++it)
134 std::cout <<
Name() <<
": Digitizing hits for " << it->first << std::endl;
137 int detID = detIDByName[it->first];
154 if(g4hits->
size() < 1)
return;
156 int detID = detIDByName[detName];
159 const PHG4Hit& g4hit = *(it->second);
162 if (! digitize_secondaries && track_id < 0)
continue;
165 double x = 0.5*(g4hit.
get_x(0) + g4hit.
get_x(1));
166 double y = 0.5*(g4hit.
get_y(0) + g4hit.
get_y(1));
167 if(!p_geomSvc->
isInPlane(detID, x, y))
continue;
169 double z = 0.5*(g4hit.
get_z(0) + g4hit.
get_z(1));
170 double px = g4hit.
get_px(0);
171 double py = g4hit.
get_py(0);
172 double pz = g4hit.
get_pz(0);
176 double x0 = x - tx*z;
177 double y0 = y - ty*z;
180 double x_ref = tx*z_ref + x0;
181 double y_ref = ty*z_ref + y0;
218 double dw = w - digiHit.
get_pos();
221 if(fabs(dw) > hodoWidth - hodoOverlap && fabs(dw) < hodoWidth)
223 if(dw < 0. && eleID != 1)
246 if(g4hits->
size() < 1)
return;
248 int detID = detIDByName[detName];
249 std::map<int, SQCalHit_v1> digiHits;
252 const PHG4Hit& g4hit = *(it->second);
255 if(digiHits.find(towerID) == digiHits.end())
277 digiHits[towerID] = digiHit;
283 for(
auto iter = digiHits.begin(); iter != digiHits.end(); ++iter)
286 iter->second.set_hit_id(digits->
size());
#define nDarkPhotonPlanes
virtual const std::string Name() const
Returns the name of this module.
@ VERBOSITY_A_LOT
Output a lot of messages.
virtual int Verbosity() const
Gets the verbosity of this module.
double getPlanePosition(int detectorID) const
static GeomSvc * instance()
singlton instance
double getPlaneOverlap(int detectorID) const
void getMeasurement(int detectorID, int elementID, double &measurement, double &dmeasurement)
Convert the detectorID and elementID to the actual hit position.
int getPlaneNElements(int detectorID) const
int getExpElementID(int detectorID, double pos_exp)
std::string getDetectorName(const int &detectorID) const
double getInterception(int detectorID, double tx, double ty, double x0, double y0) const
Get the interception of a line an a plane.
double getDCA(int detectorID, int elementID, double tx, double ty, double x0, double y0)
bool isInPlane(int detectorID, double x, double y)
See if a point is in a plane.
double getCellWidth(int detectorID) const
PHBoolean addNode(PHNode *)
Map::const_iterator ConstIterator
unsigned int size(void) const
ConstRange getHits(const unsigned int detid) const
return all hits matching a given detid
virtual float get_py(const int i) const
virtual float get_z(const int i) const
virtual float get_edep() const
virtual int get_shower_id() const
virtual float get_pz(const int i) const
virtual float get_px(const int i) const
virtual PHG4HitDefs::keytype get_hit_id() const
virtual float get_y(const int i) const
virtual float get_x(const int i) const
virtual int get_scint_id() const
virtual int get_index_l() const
virtual int get_trkid() const
PHNode * findFirst(const std::string &, const std::string &)
virtual void add_cell(short i, float edep)
virtual void set_element_id(const short id)
virtual void set_detector_id(const short a)
virtual void set_truth_y(const float a)
virtual void set_truth_pz(const float a)
virtual void set_truth_px(const float a)
virtual void set_truth_py(const float a)
virtual void set_truth_z(const float a)
virtual void set_track_id(const int a)
virtual void set_shower_id(const int a)
virtual void set_g4hit_id(const PHG4HitDefs::keytype a)
virtual void set_truth_x(const float a)
int Init(PHCompositeNode *topNode)
int InitRun(PHCompositeNode *topNode)
module initialization
SQDigitizer(const std::string &name="SQDigitizer", const int verbose=0)
void digitizePlane(const std::string &detName)
main external call, fill the digi hit vector
int process_event(PHCompositeNode *topNode)
event processing
void digitizeEMCal(const std::string &detName)
digitize the emcal hits
virtual void push_back(const SQHit *hit)=0
virtual size_t size() const =0
virtual void set_drift_distance(const float a)
virtual void set_pos(const float a)
virtual float get_pos() const
Return the absolute position of this hit. Probably the value is not properly set at present.
virtual void set_tdc_time(const float a)
virtual void set_hodo_mask(const bool a)
virtual void set_hit_id(const int a)
virtual void set_element_id(const short id)
virtual void set_detector_id(const short a)
virtual void set_in_time(const bool a)
virtual void set_pos(const float a)
virtual void set_tdc_time(const float a)
virtual void set_hodo_mask(const bool a)
virtual void set_drift_distance(const float a)
virtual void set_in_time(const bool a)
virtual void set_truth_py(const float a)
virtual void set_truth_x(const float a)
virtual void set_g4hit_id(const PHG4HitDefs::keytype a)
virtual void set_track_id(const int a)
virtual void set_truth_y(const float a)
virtual void set_truth_pz(const float a)
virtual void set_truth_z(const float a)
virtual void set_truth_px(const float a)