Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:22:04

0001 /*!
0002  *  \file       TpcPrototypeGenFitTrkFinder.C
0003  *  \brief      Refit SvtxTracks with PHGenFit.
0004  *  \details    Refit SvtxTracks with PHGenFit.
0005  *  \author     Haiwang Yu <yuhw@nmsu.edu>
0006  */
0007 
0008 #include "TpcPrototypeGenFitTrkFinder.h"
0009 
0010 #include <trackbase/TrkrCluster.h>  // for TrkrCluster
0011 #include <trackbase/TrkrClusterContainer.h>
0012 #include <trackbase/TrkrDefs.h>
0013 #include <trackreco/AssocInfoContainer.h>
0014 
0015 #include <trackbase_historic/SvtxTrackMap.h>
0016 #include <trackbase_historic/SvtxTrackMap_v1.h>
0017 #include <trackbase_historic/SvtxTrack_v1.h>
0018 #include <trackbase_historic/SvtxVertexMap.h>
0019 
0020 #include <trackbase_historic/SvtxTrack.h>
0021 #include <trackbase_historic/SvtxVertexMap_v1.h>
0022 
0023 #include <phgenfit/Fitter.h>
0024 #include <phgenfit/PlanarMeasurement.h>
0025 #include <phgenfit/Track.h>
0026 
0027 #include <fun4all/Fun4AllReturnCodes.h>
0028 #include <fun4all/PHTFileServer.h>
0029 #include <fun4all/SubsysReco.h>  // for SubsysReco
0030 
0031 #include <phool/PHCompositeNode.h>
0032 #include <phool/PHIODataNode.h>
0033 #include <phool/PHNode.h>  // for PHNode
0034 #include <phool/PHNodeIterator.h>
0035 #include <phool/PHObject.h>  // for PHObject
0036 #include <phool/getClass.h>
0037 #include <phool/phool.h>
0038 
0039 #include <phfield/PHFieldUtility.h>
0040 #include <phgeom/PHGeomUtility.h>
0041 
0042 #include <GenFit/EventDisplay.h>    // for EventDisplay
0043 #include <GenFit/RKTrackRep.h>
0044 
0045 #include <TClonesArray.h>
0046 #include <TMatrixDSymfwd.h>  // for TMatrixDSym
0047 #include <TMatrixTSym.h>     // for TMatrixTSym
0048 #include <TMatrixTUtils.h>   // for TMatrixTRow
0049 #include <TTree.h>
0050 #include <TVector3.h>
0051 
0052 #include <algorithm>                              // for find
0053 #include <set>                                    // for set
0054 
0055 #include <cassert>
0056 #include <iostream>
0057 #include <map>
0058 #include <memory>
0059 #include <utility>
0060 #include <vector>
0061 
0062 class PHField;
0063 class TGeoManager;
0064 namespace genfit
0065 {
0066 class AbsTrackRep;
0067 }
0068 namespace PHGenFit { class Measurement; }
0069 
0070 #define LogDebug(exp) std::cout << "DEBUG: " << __FILE__ << ": " << __LINE__ << ": " << exp << std::endl
0071 #define LogError(exp) std::cout << "ERROR: " << __FILE__ << ": " << __LINE__ << ": " << exp << std::endl
0072 #define LogWarning(exp) std::cout << "WARNING: " << __FILE__ << ": " << __LINE__ << ": " << exp << std::endl
0073 
0074 #define WILD_FLOAT -9999.
0075 
0076 #define _DEBUG_MODE_ 0
0077 
0078 //#define _DEBUG_
0079 
0080 using namespace std;
0081 
0082 /*
0083  * Constructor
0084  */
0085 TpcPrototypeGenFitTrkFinder::TpcPrototypeGenFitTrkFinder(const string& name, int layers)
0086   : SubsysReco(name)
0087   , _fitter(nullptr)
0088   , _track_fitting_alg_name("KalmanFitter")
0089   , nLayer(layers)
0090   , minLayer(8)
0091   , maxTracklet(50)
0092   , _primary_pid_guess(2212)
0093   , rphiWindow(2)
0094   , ZWindow(2)
0095   , _clustermap(nullptr)
0096   , _trackmap(nullptr)
0097   , _assoc_container(nullptr)
0098   , _vertexmap(nullptr)
0099   , _do_eval(false)
0100   , _eval_outname("TpcPrototypeGenFitTrkFinder_eval.root")
0101   , _eval_tree(nullptr)
0102   , _tca_ntrack(-1)
0103   , _tca_trackmap(nullptr)
0104   , _do_evt_display(false)
0105 {
0106   Verbosity(0);
0107 
0108   _event = 0;
0109 
0110   _cluster_eval_tree = nullptr;
0111   _cluster_eval_tree_x = WILD_FLOAT;
0112   _cluster_eval_tree_y = WILD_FLOAT;
0113   _cluster_eval_tree_z = WILD_FLOAT;
0114   _cluster_eval_tree_gx = WILD_FLOAT;
0115   _cluster_eval_tree_gy = WILD_FLOAT;
0116   _cluster_eval_tree_gz = WILD_FLOAT;
0117 }
0118 
0119 /*
0120  * Init
0121  */
0122 int TpcPrototypeGenFitTrkFinder::Init(PHCompositeNode* topNode)
0123 {
0124   //    CreateNodes(topNode);
0125 
0126   return Fun4AllReturnCodes::EVENT_OK;
0127 }
0128 
0129 /*
0130  * Init run
0131  */
0132 int TpcPrototypeGenFitTrkFinder::InitRun(PHCompositeNode* topNode)
0133 {
0134   CreateNodes(topNode);
0135 
0136   TGeoManager* tgeo_manager = PHGeomUtility::GetTGeoManager(topNode);
0137   PHField* field = PHFieldUtility::GetFieldMapNode(nullptr, topNode);
0138 
0139   //_fitter = new PHGenFit::Fitter("sPHENIX_Geo.root","sPHENIX.2d.root", 1.4 / 1.5);
0140   _fitter = PHGenFit::Fitter::getInstance(tgeo_manager,
0141                                           field, _track_fitting_alg_name,
0142                                           "RKTrackRep", _do_evt_display);
0143   _fitter->set_verbosity(Verbosity());
0144 
0145   if (_do_eval)
0146   {
0147     if (Verbosity() >= 1)
0148       cout << PHWHERE << " Openning file: " << _eval_outname << endl;
0149     PHTFileServer::get().open(_eval_outname, "RECREATE");
0150     init_eval_tree();
0151   }
0152 
0153   return Fun4AllReturnCodes::EVENT_OK;
0154 }
0155 /*
0156  * process_event():
0157  *  Call user instructions for every event.
0158  *  This function contains the analysis structure.
0159  *
0160  */
0161 int TpcPrototypeGenFitTrkFinder::process_event(PHCompositeNode* topNode)
0162 {
0163   _event++;
0164 
0165   if (Verbosity() > 1)
0166     std::cout << PHWHERE << "Events processed: " << _event << std::endl;
0167   //    if (_event % 1000 == 0)
0168   //        cout << PHWHERE << "Events processed: " << _event << endl;
0169 
0170   //cout << "Start TpcPrototypeGenFitTrkFinder::process_event" << endl;
0171 
0172   GetNodes(topNode);
0173 
0174   vector<std::shared_ptr<PHGenFit::Track> > rf_phgf_tracks;
0175   //  rf_phgf_tracks.clear();
0176 
0177   // Need to keep tracks if _do_evt_display
0178   if (!_do_evt_display)
0179   {
0180     rf_phgf_tracks.clear();
0181   }
0182 
0183   assert(_clustermap);
0184   set<tracklet_t> tracklets;
0185 
0186   //progressive search track finding
0187   for (int layer = 0; layer < nLayer; ++layer)
0188   {
0189     set<tracklet_t> new_tracklets(tracklets);
0190     auto cluster_range = _clustermap->getClusters(TrkrDefs::tpcId, layer);
0191 
0192     for (auto iter = cluster_range.first; iter != cluster_range.second; ++iter)
0193     {
0194       const TrkrCluster* cluster = iter->second;
0195       assert(cluster);
0196 
0197       TVector3 pos(cluster->getPosition(0), cluster->getPosition(1), cluster->getPosition(2));
0198       TVector3 n(cluster->getPosition(0), cluster->getPosition(1), 0);
0199 
0200       TVector3 n_dir(cluster->getPosition(0), cluster->getPosition(1), 0);
0201       n_dir.SetMag(1);
0202       TVector3 z_dir(0, 0, 1);
0203       TVector3 azimuth_dir(z_dir.Cross(n_dir));
0204 
0205       bool matched = false;
0206       for (const auto tracklet : tracklets)
0207       {
0208         assert(tracklet.size() >= 1);
0209         const TrkrCluster* last_cluster = tracklet.back();
0210         assert(last_cluster);
0211 
0212         TVector3 last_pos(last_cluster->getPosition(0), last_cluster->getPosition(1), last_cluster->getPosition(2));
0213 
0214         TVector3 pos_diff = pos - last_pos;
0215         const double n_residual = pos_diff.Dot(n_dir);
0216         const double z_residual = pos_diff.Dot(z_dir);
0217         const double azimuth_residual = pos_diff.Dot(azimuth_dir);
0218 
0219         assert(n_residual > 0);
0220         const double z_diff = z_residual / n_residual;
0221         const double azimuth_diff = azimuth_residual / n_residual;
0222         if (abs(azimuth_diff) < rphiWindow and abs(z_diff) < ZWindow)
0223         {  //have match
0224 
0225           if (Verbosity() >= 2)
0226           {
0227             cout << __PRETTY_FUNCTION__ << " adding cluster at layer " << layer << " to tracklet length " << tracklet.size() << endl;
0228           }
0229 
0230           auto iter_old_tracklet =
0231               new_tracklets.find(tracklet);
0232           if (iter_old_tracklet != new_tracklets.end())
0233             new_tracklets.erase(iter_old_tracklet);  //remove the shorter tracklet
0234 
0235           tracklet_t new_tracklet(tracklet);
0236           new_tracklet.push_back(cluster);
0237           new_tracklets.insert(new_tracklet);  // insert the longer tracklet
0238           matched = true;
0239 
0240           if (new_tracklets.size() >= maxTracklet)
0241           {
0242             if (Verbosity())
0243               cout << __PRETTY_FUNCTION__ << " skipping rest tracklet at layer " << layer
0244                    << " due to tracklet count " << new_tracklets.size() << " > " << maxTracklet << endl;
0245             break;
0246           }
0247         }  //have match
0248 
0249       }  //     for (auto& tracklet : tracklets)
0250 
0251       if (new_tracklets.size() >= maxTracklet)
0252       {
0253         if (Verbosity())
0254           cout << __PRETTY_FUNCTION__ << " skipping rest clusters at layer " << layer
0255                << " due to tracklet count " << new_tracklets.size() << " > " << maxTracklet << endl;
0256         break;
0257       }
0258 
0259       if (not matched)
0260       {  //new track seed from unused cluster
0261         new_tracklets.insert(tracklet_t(1, cluster));
0262         if (Verbosity() >= 2)
0263         {
0264           cout << __PRETTY_FUNCTION__ << " init tracket with cluster at layer " << layer << endl;
0265         }
0266       }  //      if (not matched)
0267 
0268     }  //      for (auto iter = cluster_range.first; iter != cluster_range.second; ++iter)
0269 
0270     tracklets = new_tracklets;
0271   }  //  for (int layer = 0; layer < nLayer; ++layer)
0272 
0273   if (Verbosity())
0274   {
0275     cout << __PRETTY_FUNCTION__ << "print initial trackets: ";
0276     for (const auto& tracklet : tracklets)
0277     {
0278       cout << ","
0279            << "size = " << tracklet.size();
0280     }
0281     cout << endl;
0282   }
0283 
0284   // track quality map
0285   multimap<double, tracklet_t> quality_tracklets_map;
0286   for (const auto& tracklet : tracklets)
0287   {
0288     double chi2Ndf = getChi2Ndf(tracklet);
0289     if (chi2Ndf > 0)
0290       quality_tracklets_map.insert(std::pair<double, tracklet_t>(chi2Ndf, tracklet));
0291   }
0292   if (Verbosity())
0293   {
0294     cout << __PRETTY_FUNCTION__ << "print fitted trackets: ";
0295     for (const auto& tracklet : quality_tracklets_map)
0296     {
0297       cout << ","
0298            << "size = " << tracklet.second.size() << " chi2/ndf = " << tracklet.first;
0299     }
0300     cout << endl;
0301   }  //  if (Verbosity())
0302 
0303   // deghost
0304   for (auto iter_current = quality_tracklets_map.begin(); iter_current != quality_tracklets_map.end(); ++iter_current)
0305   {
0306     const tracklet_t& track_current = iter_current->second;
0307 
0308     auto iter_check = iter_current;
0309     ++iter_check;
0310     for (; iter_check != quality_tracklets_map.end();)
0311     {
0312       const tracklet_t& track_check = iter_check->second;
0313 
0314       unsigned int identical_cluster = 0;
0315 
0316       for (const auto cluster : track_current)
0317       {
0318         assert(cluster);
0319 
0320         if (find(track_check.begin(), track_check.end(), cluster) != track_check.end())
0321           ++identical_cluster;
0322       }
0323 
0324       if (identical_cluster >= track_current.size() / 3 or identical_cluster >= track_check.size() / 3)
0325       {
0326         if (Verbosity())
0327         {
0328           cout << __PRETTY_FUNCTION__ << " found " << identical_cluster << "-shared ghost track size" << track_check.size() << " chi2ndf" << iter_check->first  //
0329                << " from base track size" << track_current.size() << " chi2ndf" << iter_current->first << endl;
0330         }
0331 
0332         auto iter_tmp = iter_check;
0333         ++iter_check;
0334 
0335         quality_tracklets_map.erase(iter_tmp);
0336       }
0337       else
0338       {
0339         if (Verbosity())
0340         {
0341           cout << __PRETTY_FUNCTION__ << "low " << identical_cluster << "-shared track size" << track_check.size() << " chi2ndf" << iter_check->first  //
0342                << " from base track size" << track_current.size() << " chi2ndf" << iter_current->first << endl;
0343         }
0344         ++iter_check;
0345       }
0346     }
0347   }
0348   if (Verbosity())
0349   {
0350     cout << __PRETTY_FUNCTION__ << "print deghosted trackets: ";
0351     for (const auto& tracklet : quality_tracklets_map)
0352     {
0353       cout << ","
0354            << "size = " << tracklet.second.size() << " chi2/ndf = " << tracklet.first;
0355     }
0356     cout << endl;
0357   }  //  if (Verbosity())
0358 
0359   //output
0360   assert(_trackmap);
0361   assert(_assoc_container);
0362   for (const auto& quality_tracklet : quality_tracklets_map)
0363   {
0364     const tracklet_t& tracklet = quality_tracklet.second;
0365     assert(tracklet.front());
0366     assert(tracklet.back());
0367     TVector3 pos_front(tracklet.front()->getPosition(0), tracklet.front()->getPosition(1), tracklet.front()->getPosition(2));
0368     TVector3 pos_back(tracklet.back()->getPosition(0), tracklet.back()->getPosition(1), tracklet.back()->getPosition(2));
0369 
0370     TVector3 seed_mom(pos_back - pos_front);
0371     seed_mom.SetMag(120);
0372 
0373     TVector3 seed_pos(pos_front);
0374 
0375     std::unique_ptr<SvtxTrack_v1> svtx_track(new SvtxTrack_v1());
0376 
0377     svtx_track->set_id(_trackmap->size());
0378 
0379     // dummy values, set px to make it through the minimum pT cut
0380     svtx_track->set_px(seed_mom.x());
0381     svtx_track->set_py(seed_mom.y());
0382     svtx_track->set_pz(seed_mom.z());
0383     svtx_track->set_x(seed_pos.x());
0384     svtx_track->set_y(seed_pos.y());
0385     svtx_track->set_z(seed_pos.z());
0386     for (const TrkrCluster* cluster : tracklet)
0387     {
0388       svtx_track->insert_cluster_key(cluster->getClusKey());
0389       _assoc_container->SetClusterTrackAssoc(cluster->getClusKey(), svtx_track->get_id());
0390     }
0391 
0392     _trackmap->insert(svtx_track.get());
0393   }
0394 
0395   if (_do_eval)
0396   {
0397     fill_eval_tree(topNode);
0398   }
0399   return Fun4AllReturnCodes::EVENT_OK;
0400 }
0401 
0402 double TpcPrototypeGenFitTrkFinder::getChi2Ndf(const tracklet_t& tracklet)
0403 {
0404   assert(_clustermap);
0405 
0406   if (tracklet.size() < minLayer)
0407   {
0408     if (Verbosity())
0409     {
0410       cout << __PRETTY_FUNCTION__ << " drop short tracklet size = " << tracklet.size() << " < " << minLayer << endl;
0411     }
0412     return -1;
0413   }
0414 
0415   // prepare seed
0416   assert(tracklet.front());
0417   assert(tracklet.back());
0418   TVector3 pos_front(tracklet.front()->getPosition(0), tracklet.front()->getPosition(1), tracklet.front()->getPosition(2));
0419   TVector3 pos_back(tracklet.back()->getPosition(0), tracklet.back()->getPosition(1), tracklet.back()->getPosition(2));
0420 
0421   TVector3 seed_mom(pos_back - pos_front);
0422   seed_mom.SetMag(120);
0423 
0424   TVector3 seed_pos(pos_front);
0425   TMatrixDSym seed_cov(6);
0426   for (int i = 0; i < 6; i++)
0427   {
0428     for (int j = 0; j < 6; j++)
0429     {
0430       seed_cov[i][j] = 100.;
0431     }
0432   }
0433 
0434   // Create measurements
0435   std::vector<PHGenFit::Measurement*> measurements;
0436   for (const auto cluster : tracklet)
0437   {
0438     assert(cluster);
0439     if (Verbosity() >= 2)
0440     {
0441       TrkrDefs::cluskey cluster_key = cluster->getClusKey();
0442       int layer = TrkrDefs::getLayer(cluster_key);
0443 
0444       cout << __PRETTY_FUNCTION__ << "add cluster on layer " << layer << ": ";
0445       cluster->identify();
0446     }
0447     TVector3 pos(cluster->getPosition(0), cluster->getPosition(1), cluster->getPosition(2));
0448 
0449     //TODO use u, v explicitly?
0450     TVector3 n(cluster->getPosition(0), cluster->getPosition(1), 0);
0451 
0452     PHGenFit::Measurement* meas = new PHGenFit::PlanarMeasurement(pos, n,
0453                                                                   cluster->getRPhiError(),
0454                                                                   cluster->getZError());
0455 
0456     measurements.push_back(meas);
0457 
0458   }  //    for (const auto cluster : clusterLayer)
0459 
0460   // do the fit
0461   genfit::AbsTrackRep* rep = new genfit::RKTrackRep(_primary_pid_guess);
0462   std::shared_ptr<PHGenFit::Track> phgf_track(new PHGenFit::Track(rep, seed_pos, seed_mom,
0463                                                                   seed_cov));
0464   phgf_track->addMeasurements(measurements);
0465   if (_fitter->processTrack(phgf_track.get(), false) != 0)
0466   {
0467     if (Verbosity() >= 1)
0468     {
0469       LogWarning("Track fitting failed");
0470       cout << __PRETTY_FUNCTION__ << " track->getChisq() " << phgf_track->get_chi2() << " get_ndf " << phgf_track->get_ndf()
0471            << " mom.X " << phgf_track->get_mom().X()
0472            << " mom.Y " << phgf_track->get_mom().Y()
0473            << " mom.Z " << phgf_track->get_mom().Z()
0474            << endl;
0475     }
0476     //delete track;
0477     return -1;
0478   }
0479 
0480   if (Verbosity() >= 1)
0481   {
0482     cout << __PRETTY_FUNCTION__ << " track->getChisq() " << phgf_track->get_chi2() << " get_ndf " << phgf_track->get_ndf()
0483          << " mom.X " << phgf_track->get_mom().X()
0484          << " mom.Y " << phgf_track->get_mom().Y()
0485          << " mom.Z " << phgf_track->get_mom().Z()
0486          << endl;
0487   }
0488   return phgf_track->get_chi2() / phgf_track->get_ndf();
0489 }
0490 
0491 /*
0492  * End
0493  */
0494 int TpcPrototypeGenFitTrkFinder::End(PHCompositeNode* topNode)
0495 {
0496   if (_do_eval)
0497   {
0498     if (Verbosity() >= 1)
0499       cout << PHWHERE << " Writing to file: " << _eval_outname << endl;
0500     PHTFileServer::get().cd(_eval_outname);
0501     _eval_tree->Write();
0502     _cluster_eval_tree->Write();
0503   }
0504 
0505   if (_do_evt_display)
0506   {
0507     //    if(opts[i] == 'A') drawAutoScale_ = true;
0508     //    if(opts[i] == 'B') drawBackward_ = true;
0509     //    if(opts[i] == 'D') drawDetectors_ = true;
0510     //    if(opts[i] == 'E') drawErrors_ = true;
0511     //    if(opts[i] == 'F') drawForward_ = true;
0512     //    if(opts[i] == 'H') drawHits_ = true;
0513     //    if(opts[i] == 'M') drawTrackMarkers_ = true;
0514     //    if(opts[i] == 'P') drawPlanes_ = true;
0515     //    if(opts[i] == 'S') drawScaleMan_ = true;
0516     //    if(opts[i] == 'T') drawTrack_ = true;
0517     //    if(opts[i] == 'X') drawSilent_ = true;
0518     //    if(opts[i] == 'G') drawGeometry_ = true;
0519     _fitter->getEventDisplay()->setOptions("ADEHPTG");
0520 
0521     _fitter->displayEvent();
0522   }
0523   return Fun4AllReturnCodes::EVENT_OK;
0524 }
0525 
0526 /*
0527  * dtor
0528  */
0529 TpcPrototypeGenFitTrkFinder::~TpcPrototypeGenFitTrkFinder()
0530 {
0531   delete _fitter;
0532 }
0533 
0534 /*
0535  * fill_eval_tree():
0536  */
0537 void TpcPrototypeGenFitTrkFinder::fill_eval_tree(PHCompositeNode* topNode)
0538 {
0539   //! Make sure to reset all the TTree variables before trying to set them.
0540   reset_eval_variables();
0541 
0542   int i = 0;
0543   for (SvtxTrackMap::ConstIter itr = _trackmap->begin();
0544        itr != _trackmap->end(); ++itr)
0545   {
0546     new ((*_tca_trackmap)[i])(SvtxTrack_v1)(
0547         *dynamic_cast<SvtxTrack_v1*>(itr->second));
0548   }
0549   _tca_ntrack = i;
0550 
0551   _eval_tree->Fill();
0552 
0553   return;
0554 }
0555 
0556 /*
0557  * init_eval_tree
0558  */
0559 void TpcPrototypeGenFitTrkFinder::init_eval_tree()
0560 {
0561   if (!_tca_trackmap)
0562     _tca_trackmap = new TClonesArray("SvtxTrack_v1");
0563 
0564   //! create TTree
0565   _eval_tree = new TTree("T", "TpcPrototypeGenFitTrkFinder Evaluation");
0566   _eval_tree->Branch("nTrack", &_tca_ntrack, "nTrack/I");
0567 
0568   //  _eval_tree->Branch("PrimaryParticle", _tca_particlemap);
0569   //  _eval_tree->Branch("TruthVtx", _tca_vtxmap);
0570 
0571   _eval_tree->Branch("SvtxTrack", _tca_trackmap);
0572 
0573   _cluster_eval_tree = new TTree("cluster_eval", "cluster eval tree");
0574   _cluster_eval_tree->Branch("x", &_cluster_eval_tree_x, "x/F");
0575   _cluster_eval_tree->Branch("y", &_cluster_eval_tree_y, "y/F");
0576   _cluster_eval_tree->Branch("z", &_cluster_eval_tree_z, "z/F");
0577   _cluster_eval_tree->Branch("gx", &_cluster_eval_tree_gx, "gx/F");
0578   _cluster_eval_tree->Branch("gy", &_cluster_eval_tree_gy, "gy/F");
0579   _cluster_eval_tree->Branch("gz", &_cluster_eval_tree_gz, "gz/F");
0580 }
0581 
0582 /*
0583  * reset_eval_variables():
0584  *  Reset all the tree variables to their default values.
0585  *  Needs to be called at the start of every event
0586  */
0587 void TpcPrototypeGenFitTrkFinder::reset_eval_variables()
0588 {
0589   _tca_ntrack = -1;
0590   _tca_trackmap->Clear();
0591 
0592   _cluster_eval_tree_x = WILD_FLOAT;
0593   _cluster_eval_tree_y = WILD_FLOAT;
0594   _cluster_eval_tree_z = WILD_FLOAT;
0595   _cluster_eval_tree_gx = WILD_FLOAT;
0596   _cluster_eval_tree_gy = WILD_FLOAT;
0597   _cluster_eval_tree_gz = WILD_FLOAT;
0598 }
0599 
0600 int TpcPrototypeGenFitTrkFinder::CreateNodes(PHCompositeNode* topNode)
0601 {
0602   // create nodes...
0603   PHNodeIterator iter(topNode);
0604 
0605   PHCompositeNode* dstNode = static_cast<PHCompositeNode*>(iter.findFirst(
0606       "PHCompositeNode", "DST"));
0607   if (!dstNode)
0608   {
0609     cerr << PHWHERE << "DST Node missing, doing nothing." << endl;
0610     return Fun4AllReturnCodes::ABORTEVENT;
0611   }
0612   PHNodeIterator iter_dst(dstNode);
0613 
0614   // Create the SVTX node
0615   PHCompositeNode* tb_node = dynamic_cast<PHCompositeNode*>(iter_dst.findFirst(
0616       "PHCompositeNode", "SVTX"));
0617   if (!tb_node)
0618   {
0619     tb_node = new PHCompositeNode("SVTX");
0620     dstNode->addNode(tb_node);
0621     if (Verbosity() > 0)
0622       cout << "SVTX node added" << endl;
0623   }
0624 
0625   _trackmap = new SvtxTrackMap_v1;
0626   PHIODataNode<PHObject>* tracks_node = new PHIODataNode<PHObject>(
0627       _trackmap, "SvtxTrackMap", "PHObject");
0628   tb_node->addNode(tracks_node);
0629   if (Verbosity() > 0)
0630     cout << "Svtx/SvtxTrackMap node added" << endl;
0631 
0632   _vertexmap = new SvtxVertexMap_v1;
0633   PHIODataNode<PHObject>* vertexes_node = new PHIODataNode<PHObject>(
0634       _vertexmap, "SvtxVertexMap", "PHObject");
0635   tb_node->addNode(vertexes_node);
0636   if (Verbosity() > 0)
0637     cout << "Svtx/SvtxVertexMap node added" << endl;
0638 
0639   _assoc_container = new AssocInfoContainer;
0640   PHIODataNode<PHObject>* assoc_node = new PHIODataNode<PHObject>(
0641       _assoc_container, "AssocInfoContainer", "PHObject");
0642   tb_node->addNode(assoc_node);
0643   if (Verbosity() > 0)
0644     cout << "Svtx/AssocInfoContainer node added" << endl;
0645   return Fun4AllReturnCodes::EVENT_OK;
0646 }
0647 
0648 /*
0649  * GetNodes():
0650  *  Get all the all the required nodes off the node tree
0651  */
0652 int TpcPrototypeGenFitTrkFinder::GetNodes(PHCompositeNode* topNode)
0653 {
0654   // Input Svtx Clusters
0655   //_clustermap = findNode::getClass<SvtxClusterMap>(topNode, "SvtxClusterMap");
0656   _clustermap = findNode::getClass<TrkrClusterContainer>(topNode, "TRKR_CLUSTER");
0657   if (!_clustermap && _event < 2)
0658   {
0659     return Fun4AllReturnCodes::ABORTRUN;
0660   }
0661 
0662   // Input Svtx Tracks
0663   _trackmap = findNode::getClass<SvtxTrackMap>(topNode, "SvtxTrackMap");
0664   if (!_trackmap && _event < 2)
0665   {
0666     cout << PHWHERE << " SvtxTrackMap node not found on node tree"
0667          << endl;
0668     return Fun4AllReturnCodes::ABORTEVENT;
0669   }
0670 
0671   // Input Svtx Vertices
0672   _vertexmap = findNode::getClass<SvtxVertexMap>(topNode, "SvtxVertexMap");
0673   if (!_vertexmap && _event < 2)
0674   {
0675     cout << PHWHERE << " SvtxVertexrMap node not found on node tree"
0676          << endl;
0677     return Fun4AllReturnCodes::ABORTEVENT;
0678   }
0679 
0680   return Fun4AllReturnCodes::EVENT_OK;
0681 }