File indexing completed on 2026-04-06 08:10:18
0001 #ifndef MACRO_FUN4ALLG4SPHENIX_C
0002 #define MACRO_FUN4ALLG4SPHENIX_C
0003
0004 #include "G4_Input.C"
0005 #include <G4_Global.C>
0006 #include <G4Setup_sPHENIX.C>
0007
0008 #include <Trkr_RecoInit.C>
0009 #include <Trkr_Clustering.C>
0010 #include <Trkr_TruthTables.C>
0011 #include <Trkr_Reco.C>
0012
0013 #include <phpythia8/PHPy8ParticleTrigger.h>
0014
0015 #include <decayfinder/DecayFinder.h>
0016 #include <hftrackefficiency/HFTrackEfficiency.h>
0017 #include <kfparticle_sphenix/KFParticle_sPHENIX.h>
0018
0019 #include <ffamodules/FlagHandler.h>
0020 #include <ffamodules/HeadReco.h>
0021 #include <ffamodules/SyncReco.h>
0022 #include <ffamodules/CDBInterface.h>
0023 #include <phool/PHRandomSeed.h>
0024 #include <phool/recoConsts.h>
0025
0026 #include <fun4all/Fun4AllRunNodeInputManager.h>
0027 #include <fun4all/Fun4AllDstOutputManager.h>
0028 #include <fun4all/Fun4AllOutputManager.h>
0029 #include <fun4all/Fun4AllServer.h>
0030
0031 #include <simqa_modules/QAG4SimulationTracking.h>
0032 #include <qautils/QAHistManagerDef.h>
0033
0034 #include "HF_selections.C"
0035
0036 R__LOAD_LIBRARY(libfun4all.so)
0037 R__LOAD_LIBRARY(libffamodules.so)
0038 R__LOAD_LIBRARY(libdecayfinder.so)
0039 R__LOAD_LIBRARY(libhftrackefficiency.so)
0040 R__LOAD_LIBRARY(libsimqa_modules.so)
0041
0042 int Fun4All_HFG(std::string processID = "0", std::string channel = "Kshort2pipi")
0043 {
0044 int nEvents = 4e2;
0045 std::string outDir = "./" + channel + "_20260320_DetroitMB_CR_2_mode_pTref_3p0_high_pT_boost/";
0046
0047 string makeDirectory = "mkdir -p " + outDir + "hfEff";
0048 system(makeDirectory.c_str());
0049
0050
0051 Fun4AllServer *se = Fun4AllServer::instance();
0052 se->Verbosity(1);
0053
0054 PHRandomSeed::Verbosity(1);
0055 recoConsts *rc = recoConsts::instance();
0056
0057
0058
0059 Input::PYTHIA8 = true;
0060 int particleID = 421;
0061 PYTHIA8::config_file[0] = "steeringCards/pythia8_MB_Detroit.cfg";
0062 if (channel == "Kshort2pipi")
0063 {
0064
0065
0066
0067 particleID = 310;
0068 }
0069 else if (channel == "Lambda2ppi")
0070 {
0071
0072
0073
0074 particleID = 3122;
0075 }
0076 else if (channel == "minBias")
0077 {
0078 std::cout << "Min bias simulations" << std::endl;
0079 }
0080 else
0081 {
0082 std::cout << "Your decay channel " << channel << " is not known" << std::endl;
0083 exit(1);
0084 }
0085 Input::BEAM_CONFIGURATION = Input::pp_COLLISION;
0086
0087 InputInit();
0088
0089 float abs_eta = 10;
0090
0091 PHPy8ParticleTrigger * p8_hf_signal_trigger = new PHPy8ParticleTrigger();
0092 p8_hf_signal_trigger->SetPtLow(3.);
0093 p8_hf_signal_trigger->SetPtHigh(4.);
0094 p8_hf_signal_trigger->SetEtaHighLow(abs_eta, -1*abs_eta);
0095 p8_hf_signal_trigger->SetStableParticleOnly(false);
0096
0097 if (channel != "minBias")
0098 {
0099 p8_hf_signal_trigger->AddParticles(particleID);
0100 p8_hf_signal_trigger->AddParticles(-1*particleID);
0101 p8_hf_signal_trigger->PrintConfig();
0102 INPUTGENERATOR::Pythia8[0]->register_trigger(p8_hf_signal_trigger);
0103 INPUTGENERATOR::Pythia8[0]->set_trigger_OR();
0104
0105 Input::ApplysPHENIXBeamParameter(INPUTGENERATOR::Pythia8);
0106 }
0107
0108 InputRegister();
0109
0110
0111
0112 Enable::CDB = true;
0113 rc->set_StringFlag("CDB_GLOBALTAG","ProdA_2024");
0114 rc->set_uint64Flag("TIMESTAMP",1);
0115 rc->set_IntFlag("RUNNUMBER",1);
0116
0117 Enable::MVTX_APPLYMISALIGNMENT = true;
0118 ACTSGEOM::mvtx_applymisalignment = Enable::MVTX_APPLYMISALIGNMENT;
0119
0120 SyncReco *sync = new SyncReco();
0121 se->registerSubsystem(sync);
0122
0123 HeadReco *head = new HeadReco();
0124 se->registerSubsystem(head);
0125
0126 FlagHandler *flag = new FlagHandler();
0127 se->registerSubsystem(flag);
0128
0129
0130 Enable::MBDFAKE = true;
0131 Enable::PIPE = true;
0132 Enable::PIPE_ABSORBER = true;
0133 Enable::MVTX = true;
0134 Enable::INTT = true;
0135 Enable::TPC = true;
0136 Enable::MICROMEGAS = true;
0137
0138
0139
0140 G4Init();
0141 MagnetInit();
0142 MagnetFieldInit();
0143
0144 G4Setup();
0145
0146
0147 DecayFinder *myFinder = new DecayFinder("myFinder");
0148 myFinder->Verbosity(INT_MAX);
0149 if (channel == "Kshort2pipi") myFinder->setDecayDescriptor("K_S0 -> pi^- pi^+");
0150 else myFinder->setDecayDescriptor("[Lambda0 -> proton^+ pi^-]cc");
0151 myFinder->saveDST(1);
0152 myFinder->allowPi0(1);
0153 myFinder->allowPhotons(1);
0154 myFinder->triggerOnDecay(1);
0155 myFinder->setPTmin(0.);
0156 myFinder->setEtaRange(-1*abs_eta, abs_eta);
0157 myFinder->useDecaySpecificEtaRange(false);
0158 if (channel != "minBias") se->registerSubsystem(myFinder);
0159
0160 Mbd_Reco();
0161 Mvtx_Cells();
0162 Intt_Cells();
0163 TPC_Cells();
0164 Micromegas_Cells();
0165
0166 TrackingInit();
0167
0168 Mvtx_Clustering();
0169 Intt_Clustering();
0170 TPC_Clustering();
0171 Micromegas_Clustering();
0172
0173 Tracking_Reco();
0174
0175 auto vtxfinder = new PHSimpleVertexFinder;
0176 vtxfinder->Verbosity(0);
0177 vtxfinder->setDcaCut(0.5);
0178 vtxfinder->setTrackPtCut(-99999.);
0179 vtxfinder->setBeamLineCut(1);
0180 vtxfinder->setTrackQualityCut(1000000000);
0181 vtxfinder->setNmvtxRequired(2);
0182 vtxfinder->setOutlierPairCut(0.1);
0183 se->registerSubsystem(vtxfinder);
0184
0185 Global_Reco();
0186
0187 build_truthreco_tables();
0188
0189 HFTrackEfficiency *myTrackEff = new HFTrackEfficiency("myTrackEff");
0190 myTrackEff->Verbosity(INT_MAX);
0191 myTrackEff->setDFNodeName("myFinder");
0192 myTrackEff->triggerOnDecay(1);
0193 myTrackEff->writeSelectedTrackMap(true);
0194 myTrackEff->writeOutputFile(true);
0195 std::string outputHFEffFile = outDir + "/hfEff/outputHFTrackEff_" + channel + "_" + processID + ".root";
0196 myTrackEff->setOutputFileName(outputHFEffFile);
0197 if (channel != "minBias") se->registerSubsystem(myTrackEff);
0198
0199 output_dir = outDir;
0200
0201 if (run_pipi_reco) create_hf_directories(pipi_reconstruction_name, pipi_output_dir, pipi_output_reco_file);
0202 if (run_ppi_reco) create_hf_directories(ppi_reconstruction_name, ppi_output_dir, ppi_output_reco_file);
0203
0204
0205 if (run_pipi_reco) reconstruct_pipi_mass();
0206 if (run_ppi_reco) reconstruct_ppi_mass();
0207
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217
0218
0219
0220
0221
0222
0223
0224
0225
0226
0227
0228
0229
0230
0231
0232
0233
0234
0235
0236
0237
0238
0239
0240
0241
0242
0243 se->run(nEvents);
0244
0245 se->End();
0246
0247 if (run_pipi_reco) end_kfparticle(pipi_output_reco_file, pipi_output_dir);
0248 if (run_ppi_reco) end_kfparticle(ppi_output_reco_file, ppi_output_dir);
0249
0250 gSystem->Exit(0);
0251
0252 return 0;
0253 }
0254
0255 #endif