// $Id$ // // File: DReaction_factory_Dstar.cc // Created: Thu Jun 8 10:30:22 EDT 2017 // Creator: aaustreg (on Linux halld01.jlab.org 3.10.0-514.21.1.el7.x86_64 x86_64) // #include "DReaction_factory_Dstar.h" #include "DCustomAction_CALCut.h" #include "DCustomAction_InvariantMassDifference.h" //------------------ // brun //------------------ jerror_t DReaction_factory_Dstar::brun(JEventLoop* locEventLoop, int32_t locRunNumber) { vector locBeamPeriodVector; locEventLoop->GetCalib("PHOTON_BEAM/RF/beam_period", locBeamPeriodVector); dBeamBunchPeriod = locBeamPeriodVector[0]; return NOERROR; } //------------------ // evnt //------------------ jerror_t DReaction_factory_Dstar::evnt(JEventLoop* locEventLoop, uint64_t locEventNumber) { // Make as many DReaction objects as desired DReactionStep* locReactionStep = NULL; DReaction* locReaction = NULL; //create with a unique name for each DReaction object. CANNOT (!) be "Thrown" // DOCUMENTATION: // ANALYSIS library: https://halldweb1.jlab.org/wiki/index.php/GlueX_Analysis_Software // DReaction factory: https://halldweb1.jlab.org/wiki/index.php/Analysis_DReaction /************************************************** DstarPlus Reaction Definition *************************************************/ locReaction = new DReaction("DstarPlus"); //Required: DReactionSteps to specify the channel and decay chain you want to study //Particles are of type Particle_t, an enum defined in sim-recon/src/libraries/include/particleType.h //Example: g, p -> D0, pi+, (X, p) (to look for D*+) locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(Gamma); locReactionStep->Set_TargetParticleID(Proton); locReactionStep->Add_FinalParticleID(D0); locReactionStep->Add_FinalParticleID(PiPlus); locReactionStep->Add_FinalParticleID(Unknown, true); //locReactionStep->Add_FinalParticleID(Proton, false); //true: proton missing locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak //Example: D0 -> K-, pi+ locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(D0); locReactionStep->Add_FinalParticleID(KMinus); locReactionStep->Add_FinalParticleID(PiPlus); locReactionStep->Set_KinFitConstrainInitMassFlag(true); //default: true //ignored if p4 not fit or is beam //phi, omega not constrained regardless locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak /**************************************************** DstarPlus Control Settings ****************************************************/ // Highly Recommended: Set EventStore skim query (use with "eventstore" source) // This will skip creating particle combos for events that aren't in the skims you list // Query should be comma-separated list of skims to boolean-AND together //locReaction->Set_EventStoreSkims("myskim1,myskim2,myskim3"); //boolean-AND of skims // Recommended: Type of kinematic fit to perform (default is d_NoFit) //fit types are of type DKinFitType, an enum defined in sim-recon/src/libraries/ANALYSIS/DReaction.h //Options: d_NoFit (default), d_P4Fit, d_VertexFit, d_P4AndVertexFit //P4 fits automatically constrain decaying particle masses, unless they are manually disabled locReaction->Set_KinFitType(d_P4AndVertexFit); // Highly Recommended: When generating particle combinations, reject all beam photons that match to a different RF bunch locReaction->Set_MaxPhotonRFDeltaT(0.5*dBeamBunchPeriod); //should be minimum cut value // Optional: When generating particle combinations, reject all photon candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinPhotonPIDFOM(5.73303E-7); // Optional: When generating particle combinations, reject all charged track candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinChargedPIDFOM(5.73303E-7); // Highly Recommended: Cut on number of extra "good" tracks. "Good" tracks are ones that survive the "PreSelect" (or user custom) factory. // Important: Keep cut large: Can have many ghost and accidental tracks that look "good" locReaction->Set_MaxExtraGoodTracks(4); // Highly Recommended: Enable ROOT TTree output for this DReaction // string is file name (must end in ".root"!!): doen't need to be unique, feel free to change locReaction->Enable_TTreeOutput("tree_DstarPlus.root", false); //true/false: do/don't save unused hypotheses /************************************************** DstarPlus Pre-Combo Custom Cuts *************************************************/ // Highly Recommended: Very loose invariant mass cuts, applied during DParticleComboBlueprint construction // Example: pi0 -> g, g cut // locReaction->Set_InvariantMassCut(Pi0, 0.0, 0.3); locReaction->Set_InvariantMassCut(D0, 1.80, 1.92); // Highly Recommended: Very loose DAnalysisAction cuts, applied just after creating the combination (before saving it) // Example: Missing mass of proton locReaction->Add_ComboPreSelectionAction(new DCutAction_MissingMass(locReaction, false, 2.0, 3.0)); /**************************************************** DstarPlus Analysis Actions ****************************************************/ // Recommended: Analysis actions automatically performed by the DAnalysisResults factories to histogram useful quantities. //These actions are executed sequentially, and are executed on each surviving (non-cut) particle combination //Pre-defined actions can be found in ANALYSIS/DHistogramActions_*.h and ANALYSIS/DCutActions.h //If a histogram action is repeated, it should be created with a unique name (string) to distinguish them // HISTOGRAM PID locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); // CUT PID // SYS_TOF, SYS_BCAL, SYS_FCAL, ...: DetectorSystem_t: Defined in libraries/include/GlueX.h locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, PiPlus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiPlus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiPlus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 0.75, KMinus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, KMinus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, KMinus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 3.5, Proton, SYS_FCAL)); // locReaction->Add_AnalysisAction(new DCutAction_EachPIDFOM(locReaction, 5.73303E-7)); // locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 1.0, Proton, SYS_TOF)); //cut at delta-t +/- 1.0 //false: measured data // locReaction->Add_AnalysisAction(new DCutAction_PIDTimingBeta(locReaction, 0.0, 0.9, Neutron, SYS_BCAL)); //min/max beta cut for neutrons //locReaction->Add_AnalysisAction(new DCutAction_NoPIDHit(locReaction, KPlus)); //for K+ candidates, cut tracks with no PID hit // DstarPlus deque locDstarPlus; locDstarPlus.push_back(D0); locDstarPlus.push_back(PiPlus); locReaction->Add_AnalysisAction(new DCutAction_InvariantMass(locReaction, 0, locDstarPlus, false, 1.9, 2.1)); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarPlus, false, 800, 1.9, 2.1, "DstarPlus_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PreKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarPlus, D0, false, 800, 0.1, 0.2, "PreKinFit")); // KINEMATIC FIT locReaction->Add_AnalysisAction(new DHistogramAction_KinFitResults(locReaction, 0.05)); //5% confidence level cut on pull histograms only locReaction->Add_AnalysisAction(new DCutAction_KinFitFOM(locReaction, -1)); //0% confidence level cut //require kinematic fit converges // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarPlus, false, 800, 1.9, 2.1, "DstarPlus_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PostKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarPlus, D0, false, 800, 0.1, 0.2, "PostKinFit")); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarPlus, true, 800, 1.9, 2.1, "DstarPlus_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, true, 600, 1.8, 1.92, "D0_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, true, 600, 0., 3., "KinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarPlus, D0, true, 800, 0.1, 0.2, "KinFit")); // Kinematics locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, false)); //false: measured data // locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, true, "KinFit")); //true: kinematic-fit data locReaction->Add_AnalysisAction(new DHistogramAction_TrackVertexComparison(locReaction)); locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction, "Final")); //_data.push_back(locReaction); //Register the DReaction with the factory /************************************************** DstarMinus Reaction Definition *************************************************/ locReaction = new DReaction("DstarMinus"); //Required: DReactionSteps to specify the channel and decay chain you want to study //Particles are of type Particle_t, an enum defined in sim-recon/src/libraries/include/particleType.h //Example: g, p -> D0bar, pi-, (X, p) (to look for D*+) locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(Gamma); locReactionStep->Set_TargetParticleID(Proton); locReactionStep->Add_FinalParticleID(D0); locReactionStep->Add_FinalParticleID(PiMinus); locReactionStep->Add_FinalParticleID(Unknown, true); //locReactionStep->Add_FinalParticleID(Proton, false); //true: proton missing locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak //Example: D0 -> K+, pi- locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(D0); locReactionStep->Add_FinalParticleID(KPlus); locReactionStep->Add_FinalParticleID(PiMinus); locReactionStep->Set_KinFitConstrainInitMassFlag(true); //default: true //ignored if p4 not fit or is beam //phi, omega not constrained regardless locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak /**************************************************** DstarMinus Control Settings ****************************************************/ // Highly Recommended: Set EventStore skim query (use with "eventstore" source) // This will skip creating particle combos for events that aren't in the skims you list // Query should be comma-separated list of skims to boolean-AND together //locReaction->Set_EventStoreSkims("myskim1,myskim2,myskim3"); //boolean-AND of skims // Recommended: Type of kinematic fit to perform (default is d_NoFit) //fit types are of type DKinFitType, an enum defined in sim-recon/src/libraries/ANALYSIS/DReaction.h //Options: d_NoFit (default), d_P4Fit, d_VertexFit, d_P4AndVertexFit //P4 fits automatically constrain decaying particle masses, unless they are manually disabled locReaction->Set_KinFitType(d_P4AndVertexFit); // Highly Recommended: When generating particle combinations, reject all beam photons that match to a different RF bunch locReaction->Set_MaxPhotonRFDeltaT(0.5*dBeamBunchPeriod); //should be minimum cut value // Optional: When generating particle combinations, reject all photon candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinPhotonPIDFOM(5.73303E-7); // Optional: When generating particle combinations, reject all charged track candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinChargedPIDFOM(5.73303E-7); // Highly Recommended: Cut on number of extra "good" tracks. "Good" tracks are ones that survive the "PreSelect" (or user custom) factory. // Important: Keep cut large: Can have many ghost and accidental tracks that look "good" locReaction->Set_MaxExtraGoodTracks(4); // Highly Recommended: Enable ROOT TTree output for this DReaction // string is file name (must end in ".root"!!): doen't need to be unique, feel free to change locReaction->Enable_TTreeOutput("tree_DstarMinus.root", false); //true/false: do/don't save unused hypotheses /************************************************** DstarMinus Pre-Combo Custom Cuts *************************************************/ // Highly Recommended: Very loose invariant mass cuts, applied during DParticleComboBlueprint construction // Example: pi0 -> g, g cut // locReaction->Set_InvariantMassCut(Pi0, 0.0, 0.3); locReaction->Set_InvariantMassCut(D0, 1.80, 1.92); // Highly Recommended: Very loose DAnalysisAction cuts, applied just after creating the combination (before saving it) // Example: Missing mass of proton locReaction->Add_ComboPreSelectionAction(new DCutAction_MissingMass(locReaction, false, 2.0, 3.0)); /**************************************************** DstarMinus Analysis Actions ****************************************************/ // Recommended: Analysis actions automatically performed by the DAnalysisResults factories to histogram useful quantities. //These actions are executed sequentially, and are executed on each surviving (non-cut) particle combination //Pre-defined actions can be found in ANALYSIS/DHistogramActions_*.h and ANALYSIS/DCutActions.h //If a histogram action is repeated, it should be created with a unique name (string) to distinguish them // HISTOGRAM PID locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); // CUT PID // SYS_TOF, SYS_BCAL, SYS_FCAL, ...: DetectorSystem_t: Defined in libraries/include/GlueX.h locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, PiMinus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 0.75, KPlus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, KPlus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, KPlus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_FCAL)); // locReaction->Add_AnalysisAction(new DCutAction_EachPIDFOM(locReaction, 5.73303E-7)); // locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 1.0, Proton, SYS_TOF)); //cut at delta-t +/- 1.0 //false: measured data // locReaction->Add_AnalysisAction(new DCutAction_PIDTimingBeta(locReaction, 0.0, 0.9, Neutron, SYS_BCAL)); //min/max beta cut for neutrons //locReaction->Add_AnalysisAction(new DCutAction_NoPIDHit(locReaction, KPlus)); //for K+ candidates, cut tracks with no PID hit // DstarMinus deque locDstarMinus; locDstarMinus.push_back(D0); locDstarMinus.push_back(PiMinus); locReaction->Add_AnalysisAction(new DCutAction_InvariantMass(locReaction, 0, locDstarMinus, false, 1.9, 2.1)); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarMinus, false, 800, 1.9, 2.1, "DstarMinus_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PreKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarMinus, D0, false, 800, 0.1, 0.2, "PreKinFit")); // KINEMATIC FIT locReaction->Add_AnalysisAction(new DHistogramAction_KinFitResults(locReaction, 0.05)); //5% confidence level cut on pull histograms only locReaction->Add_AnalysisAction(new DCutAction_KinFitFOM(locReaction, -1)); //0% confidence level cut //require kinematic fit converges // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarMinus, false, 800, 1.9, 2.1, "DstarMinus_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PostKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarMinus, D0, false, 800, 0.1, 0.2, "PostKinFit")); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarMinus, true, 800, 1.9, 2.1, "DstarMinus_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, true, 600, 1.8, 1.92, "D0_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, true, 600, 0., 3., "KinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarMinus, D0, true, 800, 0.1, 0.2, "KinFit")); // Kinematics locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, false)); //false: measured data // locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, true, "KinFit")); //true: kinematic-fit data locReaction->Add_AnalysisAction(new DHistogramAction_TrackVertexComparison(locReaction)); locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction, "Final")); //_data.push_back(locReaction); //Register the DReaction with the factory /************************************************** DstarZero Reaction Definition *************************************************/ locReaction = new DReaction("DstarZero"); //Required: DReactionSteps to specify the channel and decay chain you want to study //Particles are of type Particle_t, an enum defined in sim-recon/src/libraries/include/particleType.h //Example: g, p -> D0, pi0, (X, p) (to look for D*0) locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(Gamma); locReactionStep->Set_TargetParticleID(Proton); locReactionStep->Add_FinalParticleID(D0); locReactionStep->Add_FinalParticleID(Pi0); locReactionStep->Add_FinalParticleID(Unknown, true); //locReactionStep->Add_FinalParticleID(Proton, false); //true: proton missing locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak //Example: D0(bar) -> K+, pi- locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(D0); locReactionStep->Add_FinalParticleID(KPlus); locReactionStep->Add_FinalParticleID(PiMinus); locReactionStep->Set_KinFitConstrainInitMassFlag(true); //default: true //ignored if p4 not fit or is beam //phi, omega not constrained regardless locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak //Example: Pi0 -> g,g locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(Pi0); locReactionStep->Add_FinalParticleID(Gamma); locReactionStep->Add_FinalParticleID(Gamma); locReactionStep->Set_KinFitConstrainInitMassFlag(true); //default: true //ignored if p4 not fit or is beam //phi, omega not constrained regardless locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak /**************************************************** DstarZero Control Settings ****************************************************/ // Highly Recommended: Set EventStore skim query (use with "eventstore" source) // This will skip creating particle combos for events that aren't in the skims you list // Query should be comma-separated list of skims to boolean-AND together //locReaction->Set_EventStoreSkims("myskim1,myskim2,myskim3"); //boolean-AND of skims // Recommended: Type of kinematic fit to perform (default is d_NoFit) //fit types are of type DKinFitType, an enum defined in sim-recon/src/libraries/ANALYSIS/DReaction.h //Options: d_NoFit (default), d_P4Fit, d_VertexFit, d_P4AndVertexFit //P4 fits automatically constrain decaying particle masses, unless they are manually disabled locReaction->Set_KinFitType(d_P4AndVertexFit); // Highly Recommended: When generating particle combinations, reject all beam photons that match to a different RF bunch locReaction->Set_MaxPhotonRFDeltaT(0.5*dBeamBunchPeriod); //should be minimum cut value // Optional: When generating particle combinations, reject all photon candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinPhotonPIDFOM(5.73303E-7); // Optional: When generating particle combinations, reject all charged track candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinChargedPIDFOM(5.73303E-7); // Highly Recommended: Cut on number of extra "good" tracks. "Good" tracks are ones that survive the "PreSelect" (or user custom) factory. // Important: Keep cut large: Can have many ghost and accidental tracks that look "good" locReaction->Set_MaxExtraGoodTracks(4); // Highly Recommended: Enable ROOT TTree output for this DReaction // string is file name (must end in ".root"!!): doen't need to be unique, feel free to change locReaction->Enable_TTreeOutput("tree_DstarZero.root", false); //true/false: do/don't save unused hypotheses /************************************************** DstarZero Pre-Combo Custom Cuts *************************************************/ // Cut on the Minimal Photon Energy (200MeV) //locReaction->Add_ComboPreSelectionAction(new DCustomAction_CALCut(locReaction, 0.2)); // Highly Recommended: Very loose invariant mass cuts, applied during DParticleComboBlueprint construction // Example: pi0 -> g, g cut locReaction->Set_InvariantMassCut(Pi0, 0.05, 0.25); locReaction->Set_InvariantMassCut(D0, 1.80, 1.92); // Highly Recommended: Very loose DAnalysisAction cuts, applied just after creating the combination (before saving it) // Example: Missing mass of proton locReaction->Add_ComboPreSelectionAction(new DCutAction_MissingMass(locReaction, false, 2.0, 3.0)); /**************************************************** DstarZero Analysis Actions ****************************************************/ // Recommended: Analysis actions automatically performed by the DAnalysisResults factories to histogram useful quantities. //These actions are executed sequentially, and are executed on each surviving (non-cut) particle combination //Pre-defined actions can be found in ANALYSIS/DHistogramActions_*.h and ANALYSIS/DCutActions.h //If a histogram action is repeated, it should be created with a unique name (string) to distinguish them // HISTOGRAM PID locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); // CUT PID // SYS_TOF, SYS_BCAL, SYS_FCAL, ...: DetectorSystem_t: Defined in libraries/include/GlueX.h locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, PiMinus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 0.75, KPlus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, KPlus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, KPlus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 3.5, Proton, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 3.0, Gamma, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Gamma, SYS_FCAL)); // locReaction->Add_AnalysisAction(new DCutAction_EachPIDFOM(locReaction, 5.73303E-7)); // locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 1.0, Proton, SYS_TOF)); //cut at delta-t +/- 1.0 //false: measured data // locReaction->Add_AnalysisAction(new DCutAction_PIDTimingBeta(locReaction, 0.0, 0.9, Neutron, SYS_BCAL)); //min/max beta cut for neutrons //locReaction->Add_AnalysisAction(new DCutAction_NoPIDHit(locReaction, KPlus)); //for K+ candidates, cut tracks with no PID hit // DstarZero deque locDstarZero; locDstarZero.push_back(D0); locDstarZero.push_back(Pi0); locReaction->Add_AnalysisAction(new DCutAction_InvariantMass(locReaction, 0, locDstarZero, false, 1.9, 2.1)); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZero, false, 800, 1.9, 2.1, "DstarZero_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, Pi0, false, 600, 0.05, 0.25, "Pi0_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PreKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZero, D0, false, 800, 0.1, 0.2, "PreKinFit")); // KINEMATIC FIT locReaction->Add_AnalysisAction(new DHistogramAction_KinFitResults(locReaction, 0.05)); //5% confidence level cut on pull histograms only locReaction->Add_AnalysisAction(new DCutAction_KinFitFOM(locReaction, -1)); //0% confidence level cut //require kinematic fit converges // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZero, false, 800, 1.9, 2.1, "DstarZero_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, Pi0, false, 600, 0.05, 0.25, "Pi0_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PostKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZero, D0, false, 800, 0.1, 0.2, "PostKinFit")); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZero, true, 800, 1.9, 2.1, "DstarZero_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, true, 600, 1.8, 1.92, "D0_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, Pi0, true, 600, 0.05, 0.25, "Pi0_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, true, 600, 0., 3., "KinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZero, D0, true, 800, 0.1, 0.2, "KinFit")); // Kinematics locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, false)); //false: measured data // locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, true, "KinFit")); //true: kinematic-fit data locReaction->Add_AnalysisAction(new DHistogramAction_TrackVertexComparison(locReaction)); locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction, "Final")); _data.push_back(locReaction); //Register the DReaction with the factory /************************************************** DstarZeroGamma Reaction Definition *************************************************/ locReaction = new DReaction("DstarZeroGamma"); //Required: DReactionSteps to specify the channel and decay chain you want to study //Particles are of type Particle_t, an enum defined in sim-recon/src/libraries/include/particleType.h //Example: g, p -> D0, gamma, (X, p) (to look for D*0) locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(Gamma); locReactionStep->Set_TargetParticleID(Proton); locReactionStep->Add_FinalParticleID(D0); locReactionStep->Add_FinalParticleID(Gamma); locReactionStep->Add_FinalParticleID(Unknown, true); //locReactionStep->Add_FinalParticleID(Proton, false); //true: proton missing locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak //Example: D0(bar) -> K+, pi- locReactionStep = new DReactionStep(); locReactionStep->Set_InitialParticleID(D0); locReactionStep->Add_FinalParticleID(KPlus); locReactionStep->Add_FinalParticleID(PiMinus); locReactionStep->Set_KinFitConstrainInitMassFlag(true); //default: true //ignored if p4 not fit or is beam //phi, omega not constrained regardless locReaction->Add_ReactionStep(locReactionStep); dReactionStepPool.push_back(locReactionStep); //register so will be deleted later: prevent memory leak /**************************************************** DstarZeroGamma Control Settings ****************************************************/ // Highly Recommended: Set EventStore skim query (use with "eventstore" source) // This will skip creating particle combos for events that aren't in the skims you list // Query should be comma-separated list of skims to boolean-AND together //locReaction->Set_EventStoreSkims("myskim1,myskim2,myskim3"); //boolean-AND of skims // Recommended: Type of kinematic fit to perform (default is d_NoFit) //fit types are of type DKinFitType, an enum defined in sim-recon/src/libraries/ANALYSIS/DReaction.h //Options: d_NoFit (default), d_P4Fit, d_VertexFit, d_P4AndVertexFit //P4 fits automatically constrain decaying particle masses, unless they are manually disabled locReaction->Set_KinFitType(d_P4AndVertexFit); // Highly Recommended: When generating particle combinations, reject all beam photons that match to a different RF bunch locReaction->Set_MaxPhotonRFDeltaT(0.5*dBeamBunchPeriod); //should be minimum cut value // Optional: When generating particle combinations, reject all photon candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinPhotonPIDFOM(5.73303E-7); // Optional: When generating particle combinations, reject all charged track candidates with a PID confidence level < 5.73303E-7 (+/- 5-sigma) // Make sure PID errors are calculated correctly before using. //locReaction->Set_MinChargedPIDFOM(5.73303E-7); // Highly Recommended: Cut on number of extra "good" tracks. "Good" tracks are ones that survive the "PreSelect" (or user custom) factory. // Important: Keep cut large: Can have many ghost and accidental tracks that look "good" locReaction->Set_MaxExtraGoodTracks(4); // Highly Recommended: Enable ROOT TTree output for this DReaction // string is file name (must end in ".root"!!): doen't need to be unique, feel free to change locReaction->Enable_TTreeOutput("tree_DstarZeroGamma.root", false); //true/false: do/don't save unused hypotheses /************************************************** DstarZeroGamma Pre-Combo Custom Cuts *************************************************/ // Cut on the Minimal Photon Energy (200MeV) //locReaction->Add_ComboPreSelectionAction(new DCustomAction_CALCut(locReaction, 0.2)); // Highly Recommended: Very loose invariant mass cuts, applied during DParticleComboBlueprint construction // Example: pi0 -> g, g cut locReaction->Set_InvariantMassCut(D0, 1.80, 1.92); // Highly Recommended: Very loose DAnalysisAction cuts, applied just after creating the combination (before saving it) // Example: Missing mass of proton locReaction->Add_ComboPreSelectionAction(new DCutAction_MissingMass(locReaction, false, 2.0, 3.0)); /**************************************************** DstarZeroGamma Analysis Actions ****************************************************/ // Recommended: Analysis actions automatically performed by the DAnalysisResults factories to histogram useful quantities. //These actions are executed sequentially, and are executed on each surviving (non-cut) particle combination //Pre-defined actions can be found in ANALYSIS/DHistogramActions_*.h and ANALYSIS/DCutActions.h //If a histogram action is repeated, it should be created with a unique name (string) to distinguish them // HISTOGRAM PID locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction)); // CUT PID // SYS_TOF, SYS_BCAL, SYS_FCAL, ...: DetectorSystem_t: Defined in libraries/include/GlueX.h locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, PiMinus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, PiMinus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 0.75, KPlus, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, KPlus, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.0, KPlus, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_TOF)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Proton, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 3.5, Proton, SYS_FCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 3.0, Gamma, SYS_BCAL)); locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 2.5, Gamma, SYS_FCAL)); // locReaction->Add_AnalysisAction(new DCutAction_EachPIDFOM(locReaction, 5.73303E-7)); // locReaction->Add_AnalysisAction(new DCutAction_PIDDeltaT(locReaction, false, 1.0, Proton, SYS_TOF)); //cut at delta-t +/- 1.0 //false: measured data // locReaction->Add_AnalysisAction(new DCutAction_PIDTimingBeta(locReaction, 0.0, 0.9, Neutron, SYS_BCAL)); //min/max beta cut for neutrons //locReaction->Add_AnalysisAction(new DCutAction_NoPIDHit(locReaction, KPlus)); //for K+ candidates, cut tracks with no PID hit // DstarZeroGamma deque locDstarZeroGamma; locDstarZeroGamma.push_back(D0); locDstarZeroGamma.push_back(Gamma); locReaction->Add_AnalysisAction(new DCutAction_InvariantMass(locReaction, 0, locDstarZeroGamma, false, 1.9, 2.1)); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZeroGamma, false, 800, 1.9, 2.1, "DstarZero_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PreKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PreKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZeroGamma, D0, false, 800, 0.1, 0.2, "PreKinFit")); // KINEMATIC FIT locReaction->Add_AnalysisAction(new DHistogramAction_KinFitResults(locReaction, 0.05)); //5% confidence level cut on pull histograms only locReaction->Add_AnalysisAction(new DCutAction_KinFitFOM(locReaction, -1)); //0% confidence level cut //require kinematic fit converges // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZeroGamma, false, 800, 1.9, 2.1, "DstarZero_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, false, 600, 1.8, 1.92, "D0_PostKinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, false, 600, 0., 3., "PostKinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZeroGamma, D0, false, 800, 0.1, 0.2, "PostKinFit")); // HISTOGRAM MASSES //false/true: measured/kinfit data locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, 0, locDstarZeroGamma, true, 800, 1.9, 2.1, "DstarZero_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_InvariantMass(locReaction, D0, true, 600, 1.8, 1.92, "D0_KinFit")); locReaction->Add_AnalysisAction(new DHistogramAction_MissingMass(locReaction, true, 600, 0., 3., "KinFit")); locReaction->Add_AnalysisAction(new DCustomAction_InvariantMassDifference(locReaction, 0, locDstarZeroGamma, D0, true, 800, 0.1, 0.2, "KinFit")); // Kinematics locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, false)); //false: measured data // locReaction->Add_AnalysisAction(new DHistogramAction_ParticleComboKinematics(locReaction, true, "KinFit")); //true: kinematic-fit data locReaction->Add_AnalysisAction(new DHistogramAction_TrackVertexComparison(locReaction)); locReaction->Add_AnalysisAction(new DHistogramAction_PID(locReaction, "Final")); _data.push_back(locReaction); //Register the DReaction with the factory return NOERROR; } //------------------ // fini //------------------ jerror_t DReaction_factory_Dstar::fini(void) { for(size_t loc_i = 0; loc_i < dReactionStepPool.size(); ++loc_i) delete dReactionStepPool[loc_i]; //cleanup memory return NOERROR; }