Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:12:44

0001 #include "HFMLTriggerInterface.h"
0002 
0003 #include <fun4all/Fun4AllReturnCodes.h>
0004 #include <pdbcalbase/PdbParameterMap.h>
0005 #include <phool/PHCompositeNode.h>
0006 #include <phool/PHTimeServer.h>
0007 #include <phool/PHTimer.h>
0008 #include <phool/getClass.h>
0009 
0010 #include <phool/PHCompositeNode.h>
0011 
0012 #include <g4main/PHG4Hit.h>
0013 #include <g4main/PHG4Particle.h>
0014 #include <g4main/PHG4HitContainer.h>
0015 #include <g4main/PHG4TruthInfoContainer.h>
0016 
0017 //#include <trackbase/TrkrHitTruthAssoc.h>
0018 #include <trackbase/TrkrHitSetContainer.h>
0019 //#include <trackbase/TrkrClusterHitAssoc.h>
0020 #include <trackbase/TrkrHitSet.h>
0021 #include <trackbase/TrkrHit.h>
0022 //#include <trackbase/TrkrCluster.h>
0023 
0024 #include <mvtx/MvtxDefs.h>
0025 #include <mvtx/CylinderGeom_Mvtx.h>
0026 
0027 #include <g4detectors/PHG4Cell.h>
0028 #include <g4detectors/PHG4CellContainer.h>
0029 #include <g4detectors/PHG4CylinderCellGeom.h>
0030 #include <g4detectors/PHG4CylinderCellGeomContainer.h>
0031 #include <g4detectors/PHG4CylinderGeom.h>
0032 #include <g4detectors/PHG4CylinderGeomContainer.h>
0033 
0034 //#include <g4eval/SvtxEvalStack.h>
0035 
0036 #include <HepMC/GenEvent.h>
0037 #include <HepMC/GenVertex.h>
0038 #include <phhepmc/PHHepMCGenEvent.h>
0039 #include <phhepmc/PHHepMCGenEventMap.h>
0040 
0041 #include <TDatabasePDG.h>
0042 #include <TFile.h>
0043 #include <TH2F.h>
0044 #include <TH3F.h>
0045 #include <TLorentzVector.h>
0046 #include <TString.h>
0047 #include <TTree.h>
0048 
0049 #include <rapidjson/document.h>
0050 #include <rapidjson/ostreamwrapper.h>
0051 #include <rapidjson/prettywriter.h>
0052 
0053 #include <boost/format.hpp>
0054 #include <boost/property_tree/json_parser.hpp>
0055 #include <boost/property_tree/ptree.hpp>
0056 #include <boost/property_tree/xml_parser.hpp>
0057 
0058 
0059 #include <cassert>
0060 #include <cmath>
0061 #include <cstddef>
0062 #include <iostream>
0063 
0064 using namespace std;
0065 
0066 HFMLTriggerInterface::HFMLTriggerInterface(std::string filename)
0067   : SubsysReco("HFMLTriggerInterface")
0068   , _ievent(0)
0069   , _f(nullptr)
0070   , _eta_min(-1)
0071   , _eta_max(1)
0072   , _embedding_id(1)
0073   , _nlayers_maps(3)
0074   , m_hitsets(nullptr)
0075   , m_GenEventMap(nullptr)
0076   , m_truthInfo(nullptr)
0077   , m_g4hits_mvtx(nullptr)
0078  // , _svtxevalstack(nullptr)
0079  // , m_hit_truth_map(nullptr)
0080   , m_Flags(nullptr)
0081   , m_Geoms(nullptr)
0082   , m_hitStaveLayer(nullptr)
0083   , m_hitModuleHalfStave(nullptr)
0084   , m_hitChipModule(nullptr)
0085   , m_hitLayerMap(nullptr)
0086   , m_hitPixelPhiMap(nullptr)
0087   , m_hitPixelPhiMapHL(nullptr)
0088   , m_hitPixelZMap(nullptr)
0089 {
0090   _foutname = filename;
0091 }
0092 
0093 int HFMLTriggerInterface::Init(PHCompositeNode* topNode)
0094 {
0095   _ievent = 0;
0096 
0097   _f = new TFile((_foutname + string(".root")).c_str(), "RECREATE");
0098 
0099   m_jsonOut.open((_foutname + string(".json")).c_str(), fstream::out);
0100 
0101   m_jsonOut << "{" << endl;
0102   m_jsonOut << "\"Events\" : [" << endl;
0103 
0104   //  _h2 = new TH2D("h2", "", 100, 0, 100.0, 40, -2, +2);
0105   //  _h2_b = new TH2D("h2_b", "", 100, 0, 100.0, 40, -2, +2);
0106   //  _h2_c = new TH2D("h2_c", "", 100, 0, 100.0, 40, -2, +2);
0107   //  _h2all = new TH2D("h2all", "", 100, 0, 100.0, 40, -2, +2);
0108 
0109   m_hitLayerMap = new TH3F("hitLayerMap", "hitLayerMap", 600, -10, 10, 600, -10, 10, 10, -.5, 9.5);
0110   m_hitLayerMap->SetTitle("hitLayerMap;x [mm];y [mm];Half Layers");
0111 
0112   m_hitPixelPhiMap = new TH3F("hitPixelPhiMap", "hitPixelPhiMap", 16000, -.5, 16000 - .5, 600, -M_PI, M_PI, 10, -.5, 9.5);
0113   m_hitPixelPhiMap->SetTitle("hitPixelPhiMap;PixelPhiIndex in layer;phi [rad];Half Layers Index");
0114   m_hitPixelPhiMapHL = new TH3F("hitPixelPhiMapHL", "hitPixelPhiMapHL", 16000, -.5, 16000 - .5, 600, -M_PI, M_PI, 10, -.5, 9.5);
0115   m_hitPixelPhiMapHL->SetTitle("hitPixelPhiMap;PixelPhiIndex in half layer;phi [rad];Half Layers Index");
0116   m_hitPixelZMap = new TH3F("hitPixelZMap", "hitPixelZMap", 16000, -.5, 16000 - .5, 600, 15, 15, 10, -.5, 9.5);
0117   m_hitPixelZMap->SetTitle("hitPixelZMap;hitPixelZMap;z [cm];Half Layers");
0118 
0119   m_hitStaveLayer = new TH2F("hitStaveLayer", "hitStaveLayer", 100, -.5, 100 - .5, 10, -.5, 9.5);
0120   m_hitStaveLayer->SetTitle("hitStaveLayer;Stave index;Half Layers");
0121   m_hitModuleHalfStave = new TH2F("hitModuleHalfStave", "hitModuleHalfStave", 100, -.5, 100 - .5, 10, -.5, 9.5);
0122   m_hitModuleHalfStave->SetTitle("hitModuleHalfStave;Module index;Half Stave");
0123   m_hitChipModule = new TH2F("hitChipModule", "hitChipModule", 100, -.5, 100 - .5, 10, -.5, 9.5);
0124   m_hitChipModule->SetTitle("hitChipModule;Chip;Module");
0125 
0126   return Fun4AllReturnCodes::EVENT_OK;
0127 }
0128 
0129 int HFMLTriggerInterface::InitRun(PHCompositeNode* topNode)
0130 {
0131   m_Geoms =
0132       findNode::getClass<PHG4CylinderGeomContainer>(topNode, "CYLINDERGEOM_MVTX");
0133   if (!m_Geoms)
0134   {
0135     std::cout << PHWHERE << "ERROR: Can't find node CYLINDERCELLGEOM_MVTX" << std::endl;
0136     return Fun4AllReturnCodes::ABORTRUN;
0137   }
0138 
0139   m_truthInfo = findNode::getClass<PHG4TruthInfoContainer>(topNode, "G4TruthInfo");
0140   if (!m_truthInfo)
0141   {
0142     std::cout << PHWHERE << "ERROR: Can't find node G4TruthInfo" << std::endl;
0143     return Fun4AllReturnCodes::ABORTRUN;
0144   }
0145 
0146   m_Flags = findNode::getClass<PdbParameterMap>(topNode, "HFMLTrigger_HepMCTriggerFlags");
0147   if (!m_Flags)
0148   {
0149     cout << "HFMLTriggerInterface::InitRun - WARNING - missing HFMLTrigger_HepMCTriggerFlags" << endl;
0150   }
0151 
0152   return Fun4AllReturnCodes::EVENT_OK;
0153 }
0154 
0155 int HFMLTriggerInterface::process_event(PHCompositeNode* topNode)
0156 {
0157   // load nodes
0158   auto res = load_nodes(topNode);
0159   if (res != Fun4AllReturnCodes::EVENT_OK)
0160     return res;
0161 
0162 
0163  // if (!_svtxevalstack)
0164  // {
0165  //   _svtxevalstack = new SvtxEvalStack(topNode);
0166  //   _svtxevalstack->set_strict(0);
0167  //   _svtxevalstack->set_verbosity(Verbosity() + 1);
0168  // }
0169  // else
0170  // {
0171  //   _svtxevalstack->next_event(topNode);
0172  // }
0173 
0174   //  SvtxVertexEval* vertexeval = _svtxevalstack->get_vertex_eval();
0175   //  SvtxTrackEval* trackeval = _svtxevalstack->get_track_eval();
0176  // SvtxClusterEval* clustereval = _svtxevalstack->get_cluster_eval();
0177   //SvtxHitEval* hiteval = _svtxevalstack->get_hit_eval();
0178   //    SvtxTruthEval* trutheval = _svtxevalstack->get_truth_eval();
0179 
0180   assert(m_GenEventMap);
0181 
0182   PHHepMCGenEvent* genevt = m_GenEventMap->get(_embedding_id);
0183   if (!genevt)
0184   {
0185     std::cout << PHWHERE << " - Fatal error - node PHHepMCGenEventMap missing subevent with embedding ID " << _embedding_id;
0186     std::cout << ". Print PHHepMCGenEventMap:";
0187     m_GenEventMap->identify();
0188     return Fun4AllReturnCodes::ABORTEVENT;
0189   }
0190 
0191   HepMC::GenEvent* theEvent = genevt->getEvent();
0192   assert(theEvent);
0193   if (Verbosity())
0194   {
0195     cout << "HFMLTriggerInterface::process_event - process HepMC::GenEvent with signal_process_id = "
0196          << theEvent->signal_process_id();
0197     if (theEvent->signal_process_vertex())
0198     {
0199       cout << " and signal_process_vertex : ";
0200       theEvent->signal_process_vertex()->print();
0201     }
0202     cout << "  - Event record:" << endl;
0203     theEvent->print();
0204   }
0205 
0206   // property tree preparation
0207   //  using boost::property_tree::ptree;
0208   //  ptree pt;
0209 
0210   rapidjson::Document d;
0211   d.SetObject();
0212   rapidjson::Document::AllocatorType& alloc = d.GetAllocator();
0213 
0214   auto loadCoordinate = [&](double x, double y, double z) {
0215     //    ptree vertexTree;
0216     rapidjson::Value vertexTree(rapidjson::kArrayType);
0217 
0218     //    ptree vertexX;
0219     //    vertexX.put("", x);
0220     //    vertexTree.push_back(make_pair("", vertexX));
0221 
0222     //    ptree vertexY;
0223     //    vertexY.put("", y);
0224     //    vertexTree.push_back(make_pair("", vertexY));
0225 
0226     //    ptree vertexZ;
0227     //    vertexZ.put("", z);
0228     //    vertexTree.push_back(make_pair("", vertexZ));
0229 
0230     vertexTree.PushBack(x, alloc).PushBack(y, alloc).PushBack(z, alloc);
0231 
0232     return vertexTree;
0233   };
0234 
0235   // Create a root
0236   //  ptree pTree;
0237 
0238   // meta data
0239   //  ptree metaTree;
0240   rapidjson::Value metaTree(rapidjson::kObjectType);
0241 
0242   metaTree.AddMember("Description", "These are meta data for this event. Not intended to use in ML algorithm", alloc);
0243   metaTree.AddMember("EventID", _ievent, alloc);
0244   metaTree.AddMember("Unit", "cm", alloc);
0245   metaTree.AddMember("CollisionVertex",
0246                      loadCoordinate(genevt->get_collision_vertex().x(),
0247                                     genevt->get_collision_vertex().y(),
0248                                     genevt->get_collision_vertex().z()),
0249                      alloc);
0250 
0251   metaTree.AddMember("Layer_Count", _nlayers_maps, alloc);
0252   metaTree.AddMember("PixelHalfLayerIndex_Count", _nlayers_maps * 2, alloc);
0253 
0254   for (unsigned int layer = 0; layer < _nlayers_maps; ++layer)
0255   {
0256     CylinderGeom_Mvtx* geom = dynamic_cast<CylinderGeom_Mvtx*>(m_Geoms->GetLayerGeom(layer));
0257     assert(geom);
0258 
0259     //    ptree layerDescTree;
0260     rapidjson::Value layerDescTree(rapidjson::kObjectType);
0261 
0262     static const unsigned int nChip(9);
0263 
0264     layerDescTree.AddMember("PixelPhiIndexInLayer_Count", geom->get_N_staves() * geom->get_NX(), alloc);
0265     layerDescTree.AddMember("PixelPhiIndexInHalfLayer_Count", geom->get_N_staves() * geom->get_NX() / 2, alloc);
0266     layerDescTree.AddMember("PixelZIndex_Count", nChip * geom->get_NZ(), alloc);
0267     layerDescTree.AddMember("HalfLayer_Count", 2, alloc);
0268     layerDescTree.AddMember("Stave_Count", geom->get_N_staves(), alloc);
0269     layerDescTree.AddMember("Chip_Count", nChip, alloc);
0270     layerDescTree.AddMember("Pixel_Count", geom->get_NX() * geom->get_NZ(), alloc);
0271 
0272     //    metaTree.AddMember(
0273     //        str(boost::format{"Layer%1%"} % layer).c_str(),
0274     //        layerDescTree, alloc);
0275     rapidjson::Value keyName(str(boost::format{"Layer%1%"} % layer).c_str(), alloc);
0276     metaTree.AddMember(keyName,
0277                        layerDescTree, alloc);
0278   }
0279 
0280   //  ptree truthTriggerFlagTree;
0281   rapidjson::Value truthTriggerFlagTree(rapidjson::kObjectType);
0282 
0283   truthTriggerFlagTree.AddMember("Description",
0284                                  "These are categorical true/false MonteCalo truth tags for the event. These are only known in training sample. This would be trigger output in real data processing.",
0285                                  alloc);
0286   //    truthTriggerFlagTree.AddMember("ExampleSignal1", true, alloc);
0287   //    truthTriggerFlagTree.AddMember("ExampleSignal2", false, alloc);
0288   rapidjson::Value flagsTree(rapidjson::kObjectType);
0289   if (m_Flags)
0290   {
0291     auto range = m_Flags->get_iparam_iters();
0292 
0293     for (auto flagIter = range.first; flagIter != range.second; ++flagIter)
0294     {
0295 //      rapidjson::Value aFlag(rapidjson::kObjectType);
0296 
0297       const string& name = flagIter->first;
0298       rapidjson::Value keyName(name.c_str(), alloc);
0299       const bool flag = flagIter->second > 0 ? true : false;
0300 
0301       flagsTree.AddMember(keyName, flag, alloc);
0302 //      flagsTree.PushBack(aFlag, alloc);
0303     }
0304   }
0305   truthTriggerFlagTree.AddMember("Flags", flagsTree, alloc);
0306 
0307   // Raw hits
0308   //  ptree rawHitTree;
0309   rapidjson::Value rawHitTree(rapidjson::kObjectType);
0310   rawHitTree.AddMember("Description",
0311                        "Raw data in the event in an unordered set of hit ID. To help in visualization stage, the coordinate of the hit is also appended. These would be input in real data.",
0312                        alloc);
0313 
0314   //  rawHitTree.put("LayerRage", "0-2");
0315 
0316   //  ptree rawHitsTree;
0317   rapidjson::Value rawHitsTree(rapidjson::kArrayType);
0318 
0319   assert(m_hitsets);
0320   unsigned int hitID(0);
0321   auto hitset_range = m_hitsets->getHitSets(TrkrDefs::TrkrId::mvtxId);
0322   for (auto hitset_iter =  hitset_range.first; hitset_iter != hitset_range.second; ++hitset_iter)
0323   {
0324     auto hitsetkey = hitset_iter->first;
0325     auto hit_range = hitset_iter->second->getHits();
0326     for (auto hit_iter = hit_range.first;
0327          hit_iter != hit_range.second;
0328          ++hit_iter)
0329     {
0330       auto hitkey = hit_iter->first;
0331       auto hit = hit_iter->second;
0332       assert(hit);
0333       unsigned int layer = TrkrDefs::getLayer(hitset_iter->first);
0334       if (layer < _nlayers_maps)
0335       {
0336         CylinderGeom_Mvtx* geom = dynamic_cast<CylinderGeom_Mvtx*>(m_Geoms->GetLayerGeom(layer));
0337         assert(geom);
0338 
0339         unsigned int pixel_x = MvtxDefs::getRow(hitkey);
0340         unsigned int pixel_z = MvtxDefs::getCol(hitkey);
0341         assert((int)pixel_x < geom->get_NX());
0342         assert((int)pixel_z < geom->get_NZ());
0343         unsigned int chip = MvtxDefs::getChipId(hitsetkey);
0344         unsigned int stave = MvtxDefs::getStaveId(hitsetkey);
0345         TVector3 local_coords = geom->get_local_coords_from_pixel(pixel_x, pixel_z);
0346         TVector3 world_coords = geom->get_world_from_local_coords(stave,
0347                                                                   0,
0348                                                                   0,
0349                                                                   chip,
0350                                                                   local_coords);
0351 
0352         //unsigned int halflayer = (int) pixel_x >= geom->get_NX() / 2 ? 0 : 1;
0353 
0354         //unsigned int halfLayerIndex(layer * 2 + halflayer);
0355         //unsigned int pixelPhiIndex(
0356         //  cell->get_stave_index() * geom->get_NX() + pixel_x);
0357         //unsigned int pixelPhiIndexHL(
0358         //  cell->get_stave_index() * geom->get_NX() / 2 + pixel_x % (geom->get_NX() / 2));
0359         // unsigned int pixelZIndex(cell->get_chip_index() * geom->get_NZ() + pixel_z);
0360 
0361         //      ptree hitTree;
0362         rapidjson::Value hitTree(rapidjson::kObjectType);
0363 
0364         //      ptree hitIDTree;
0365         rapidjson::Value hitIDTree(rapidjson::kObjectType);
0366         hitIDTree.AddMember("HitSequenceInEvent", hitID, alloc);
0367 
0368         //hitIDTree.AddMember("PixelHalfLayerIndex", halfLayerIndex, alloc);
0369         //hitIDTree.AddMember("PixelPhiIndexInLayer", pixelPhiIndex, alloc);
0370         //hitIDTree.AddMember("PixelPhiIndexInHalfLayer", pixelPhiIndexHL, alloc);
0371         //hitIDTree.AddMember("PixelZIndex", pixelZIndex, alloc);
0372 
0373         hitIDTree.AddMember("Layer", layer, alloc);
0374         //hitIDTree.AddMember("HalfLayer", halflayer, alloc);
0375         hitIDTree.AddMember("Stave", stave, alloc);
0376         //      hitIDTree.put("HalfStave", cell->get_half_stave_index());
0377         //      hitIDTree.put("Module", cell->get_module_index());
0378         hitIDTree.AddMember("Chip", chip, alloc);
0379         hitIDTree.AddMember("Pixel_x", pixel_x, alloc);
0380         hitIDTree.AddMember("Pixel_z", pixel_z, alloc);
0381         hitTree.AddMember("ID", hitIDTree, alloc);
0382 
0383         hitTree.AddMember("Coordinate",
0384                           loadCoordinate(world_coords.x(),
0385                                          world_coords.y(),
0386                                          world_coords.z()),
0387                           alloc);
0388 
0389         //      rawHitsTree.add_child("MVTXHit", hitTree);
0390         rawHitsTree.PushBack(hitTree, alloc);
0391 
0392         m_hitStaveLayer->Fill(stave, layer);
0393         m_hitModuleHalfStave->Fill(stave, layer);
0394         m_hitChipModule->Fill(chip, stave);
0395 
0396         m_hitLayerMap->Fill(world_coords.x(), world_coords.y(), layer);
0397         //m_hitPixelPhiMap->Fill(pixelPhiIndex, atan2(world_coords.y(), world_coords.x()), layer);
0398         //m_hitPixelPhiMapHL->Fill(pixelPhiIndexHL, atan2(world_coords.y(), world_coords.x()), layer);
0399         //m_hitPixelZMap->Fill(pixelZIndex, world_coords.z(), halfLayerIndex);
0400         ++hitID;
0401       }  //    if (layer < _nlayers_maps)
0402     }  //   for (TrkrHitSetContainer::ConstIterator hitset_iter = hitset_range.first;
0403   }  //   for (TrkrHitSet::ConstIterator hit_iter = hit_range.first;
0404   rawHitTree.AddMember("MVTXHits", rawHitsTree, alloc);
0405 
0406   // Truth hits
0407   //  ptree truthHitTree;
0408   rapidjson::Value truthHitTree(rapidjson::kObjectType);
0409   truthHitTree.AddMember("Description", "From the MonteCalo truth information, pairs of track ID and subset of RawHit that belong to the track. These are not presented in real data. The track ID is arbitary.",
0410                          alloc);
0411 
0412     assert(m_g4hits_mvtx);
0413     std::multimap<const int, const unsigned long long> m_track_g4hits_map;
0414     for (auto g4hit_iter = m_g4hits_mvtx->getHits().first;
0415         g4hit_iter != m_g4hits_mvtx->getHits().second;
0416         ++g4hit_iter)
0417     {
0418       PHG4Hit *g4hit = g4hit_iter->second;
0419       m_track_g4hits_map.insert(make_pair(g4hit->get_trkid(), g4hit_iter->first));
0420     }
0421 
0422   //  ptree truthTracksTree;
0423   rapidjson::Value truthTracksTree(rapidjson::kArrayType);
0424 
0425   // get set of primary particle
0426   assert(m_truthInfo);
0427   auto pp_range = m_truthInfo->GetPrimaryParticleRange();
0428   for (auto pp_iter = pp_range.first;
0429        pp_iter != pp_range.second;
0430        ++pp_iter)
0431   {
0432     PHG4Particle *g4particle = pp_iter->second;
0433     assert(g4particle);
0434 
0435     if (! m_track_g4hits_map.count(g4particle->get_track_id()))
0436     {
0437       if (Verbosity() >= VERBOSITY_MORE)
0438         std::cout << "WARNING: G4 particle " << g4particle->get_track_id() \
0439                   << " does not hit any MVTX layer." << std::endl;
0440       continue;
0441     }
0442 
0443     //    ptree trackHitTree;
0444     rapidjson::Value trackHitTree(rapidjson::kArrayType);
0445     //unsigned int nMAPS(0);
0446     auto g4hits_iter = m_track_g4hits_map.equal_range(g4particle->get_track_id());
0447     for (auto& g4hit_iter = g4hits_iter.first;
0448         g4hit_iter != g4hits_iter.second; ++g4hit_iter)
0449     {
0450       //          ptree hitIDTree;
0451       //          hitIDTree.put("", g4hit_key);
0452       trackHitTree.PushBack(static_cast<uint64_t>(g4hit_iter->second), alloc);
0453     }  //   for (auto& clus_key : clus_keys)
0454 
0455     //      ptree trackTree;
0456     rapidjson::Value trackTree(rapidjson::kObjectType);
0457     trackTree.AddMember("TrackSequenceInEvent", g4particle->get_track_id(), alloc);
0458     trackTree.AddMember("HitSequenceInEvent", trackHitTree, alloc);
0459 
0460     trackTree.AddMember("ParticleTypeID", g4particle->get_pid(), alloc);
0461     trackTree.AddMember("TrackMomentum",
0462                         loadCoordinate(g4particle->get_px(),
0463                                        g4particle->get_py(),
0464                                        g4particle->get_pz()),
0465                         alloc);
0466       //      trackTree.put("TrackDCA3DXY", track->get_dca3d_xy());
0467       //      trackTree.put("TrackDCA3DZ", track->get_dca3d_z());
0468 
0469       //      trackTree.add_child("TruthHit", trackHitTree);
0470 
0471       //      truthTracksTree.add_child("TruthTrack", trackTree);
0472     truthTracksTree.PushBack(trackTree, alloc);
0473   }  //  for (PHG4TruthInfoContainer::ConstIterator iter = range.first;
0474 
0475   truthHitTree.AddMember("TruthTracks", truthTracksTree, alloc);
0476 
0477   //output
0478   d.AddMember("MetaData", metaTree, alloc);
0479   d.AddMember("TruthTriggerFlag", truthTriggerFlagTree, alloc);
0480   d.AddMember("RawHit", rawHitTree, alloc);
0481   d.AddMember("TruthHit", truthHitTree, alloc);
0482 
0483   assert(m_jsonOut.is_open());
0484 
0485   if (_ievent > 0)
0486   {
0487     m_jsonOut << "," << endl;
0488   }
0489 
0490   //  write_json(m_jsonOut, pTree);
0491   //  write_xml(m_jsonOut, jsonTree);
0492 
0493   //  d.AddMember("Test", 1, d.GetAllocator());
0494 
0495   rapidjson::OStreamWrapper osw(m_jsonOut);
0496   rapidjson::PrettyWriter<rapidjson::OStreamWrapper> writer(osw);
0497 
0498   d.Accept(writer);
0499 
0500   ++_ievent;
0501 
0502   return Fun4AllReturnCodes::EVENT_OK;
0503 }
0504 
0505 int HFMLTriggerInterface::End(PHCompositeNode* topNode)
0506 {
0507   if (_f)
0508   {
0509     _f->cd();
0510     _f->Write();
0511   }
0512 
0513   if (m_jsonOut.is_open())
0514   {
0515     m_jsonOut << "]" << endl;
0516     m_jsonOut << "}" << endl;
0517 
0518     m_jsonOut.close();
0519   }
0520 
0521   cout << "HFMLTriggerInterface::End - output to " << _foutname << ".*" << endl;
0522 
0523   return Fun4AllReturnCodes::EVENT_OK;
0524 }
0525 
0526 int HFMLTriggerInterface::load_nodes(PHCompositeNode* topNode)
0527 {
0528 
0529   m_hitsets = findNode::getClass<TrkrHitSetContainer>(topNode, "TRKR_HITSET");
0530   if (!m_hitsets)
0531   {
0532     std::cout << PHWHERE << "ERROR: Can't find node TRKR_HITSET" << std::endl;
0533     return Fun4AllReturnCodes::ABORTEVENT;
0534   }
0535 
0536   //m_hit_truth_map = findNode::getClass<TrkrHitTruthAssoc>(topNode, "TRKR_HITTRUTHASSOC");
0537   //if (!m_hit_truth_map)
0538   //{
0539   //  std::cout << PHWHERE << " unable to find DST node TRKR_HITTRUTHASSOC" << std::endl;
0540   //  return Fun4AllReturnCodes::ABORTEVENT;
0541   //}
0542 
0543   //m_cluster_hit_map = findNode::getClass<TrkrClusterHitAssoc>(topNode, "TRKR_CLUSTERHITASSOC");
0544   //if (!m_cluster_hit_map)
0545   //{
0546   //  std::cout << PHWHERE << " unable to find DST node TRKR_CLUSTERHITASSOC" << std::endl;
0547   //  return Fun4AllReturnCodes::ABORTEVENT;
0548   //}
0549 
0550   m_g4hits_mvtx = findNode::getClass<PHG4HitContainer>(topNode, "G4HIT_MVTX");
0551   if (!m_g4hits_mvtx)
0552   {
0553     std::cout << PHWHERE << " unable to find DST node G4HIT_MVTX" << std::endl;
0554     return Fun4AllReturnCodes::ABORTEVENT;
0555   }
0556 
0557   m_GenEventMap = findNode::getClass<PHHepMCGenEventMap>(topNode, "PHHepMCGenEventMap");
0558   if (!m_GenEventMap)
0559   {
0560     std::cout << PHWHERE << " - Fatal error - missing node PHHepMCGenEventMap" << std::endl;
0561     return Fun4AllReturnCodes::ABORTEVENT;
0562   }
0563 
0564   return Fun4AllReturnCodes::EVENT_OK;
0565 }