GENIEGenerator
Loading...
Searching...
No Matches
gCalcHEDISDiffXsec.cxx
Go to the documentation of this file.
1//____________________________________________________________________________
2/*!
3\program gcalchedisdiffxsec
4\brief GENIE utility program calculating differential cross sections from
5 HEDIS package.
6 Syntax :
7 gcalchedisdiffxsec -p nu -t tgt -o root_file
8 -x table_type
9 --tune genie_tune
10 --event-generator-list list_name
11 Note :
12 [] marks optional arguments.
13 <> marks a list of arguments out of which only one can be
14 selected at any given time.
15 Options :
16 -p
17 the neutrino pdg code
18 -t
19 the target pdg code (format: 10LZZZAAAI)
20 -x
21 path to ascii file with x values
22 -y
23 path to ascii file with y values
24 -e
25 path to ascii file with energy values
26 -o
27 output ROOT file name
28 --tune
29 Specifies a GENIE comprehensive neutrino interaction model tune.
30 --event-generator-list
31 List of event generators to load in event generation drivers.
32 [default: "Default"].
33 --message-thresholds
34 Allows users to customize the message stream thresholds.
35 The thresholds are specified using an XML file.
36 See $GENIE/config/Messenger.xml for the XML schema.
37 *** See the User Manual for more details and examples. ***
38\author Alfonso Garcia <alfonsog \at nikhef.nl>
39 NIKHEF (Amsterdam)
40\created May 26, 2021
41\cpright Copyright (c) 2003-2025, The GENIE Collaboration
42 For the full text of the license visit http://copyright.genie-mc.org
43 or see $GENIE/LICENSE
44*/
45//____________________________________________________________________________
46
47#include "TMath.h"
48#include "TSystem.h"
49#include "TTree.h"
50#include "TFile.h"
51
65
66#include <fstream>
67#include <iostream>
68
69using namespace genie;
70
71void PrintSyntax (void);
72void DecodeCommandLine (int argc, char * argv[]);
73
74vector<double> ReadListFromPath (string path);
75
76const double epsilon = 1e-5;
77
78//----------------------------------------------------------------------------------------
79//----------------------------------------------------------------------------------------
80//INPUT ARGUMENTS
81//----------------------------------------------------------------------------------------
82//----------------------------------------------------------------------------------------
83
84int fNu = -1;
85int fTgt = -1;
86string fPathXlist = "";
87string fPathYlist = "";
88string fPathElist = "";
89string fOutFileName = "";
90bool fSaveAll = false;
91
92//----------------------------------------------------------------------------------------
93//----------------------------------------------------------------------------------------
94//MAIN PROGRAM
95//----------------------------------------------------------------------------------------
96//----------------------------------------------------------------------------------------
97
98int main(int argc, char ** argv) {
99
100 DecodeCommandLine(argc,argv);
101
103 utils::app_init::MesgThresholds(RunOpt::Instance()->MesgThresholdFiles());
104
105 string fChannel = RunOpt::Instance()->EventGeneratorList();
106
107 GEVGDriver evg_driver;
108 InitialState init_state(fTgt, fNu);
109 evg_driver.SetEventGeneratorList(fChannel);
110 evg_driver.Configure(init_state);
111
112 vector<double> ve = ReadListFromPath(fPathElist);
113 vector<double> vx = ReadListFromPath(fPathXlist);
114 vector<double> vy = ReadListFromPath(fPathYlist);
115
116 double widthx = (TMath::Log10(vx[1])-TMath::Log10(vx[0]));
117
118 LOG("gcalchedisdiffxsec", pINFO) << widthx;
119
120 int quark;
121 double ei;
122 double bx;
123 double by;
124 double dxsec;
125
126 TString treename = Form("diffxsec_nu%d_tgt%d_%s",fNu,fTgt,fChannel.c_str());
127
128 LOG("gcalchedisdiffxsec", pINFO) << treename;
129
130 TTree * tree = new TTree(treename,treename);
131 tree->Branch( "Quark", &quark, "Quark/I" );
132 tree->Branch( "Ei", &ei, "Ei/D" );
133 tree->Branch( "Bx", &bx, "Bx/D" );
134 tree->Branch( "By", &by, "By/D" );
135 tree->Branch( "DiffXsec", &dxsec, "DiffXsec/D" );
136
137 const InteractionList * intlst = evg_driver.Interactions();
138
139 InteractionList::const_iterator intliter;
140 for(intliter = intlst->begin(); intliter != intlst->end(); ++intliter) {
141
142 const Interaction * interaction = *intliter;
143
144 quark = 10000 * interaction->InitState().Tgt().HitQrkPdg();
145 if (!interaction->InitState().Tgt().HitSeaQrk()) {
146 if ( quark>0 ) quark += 100;
147 else quark -= 100;
148 }
149 quark += 1 * interaction->ExclTag().FinalQuarkPdg();
150
151 LOG("gcalchedisdiffxsec", pINFO) << "Current interaction: " << interaction->AsString();
152
153 const XSecAlgorithmI * xsec_alg = evg_driver.FindGenerator(interaction)->CrossSectionAlg();
154
155 for ( unsigned i=0; i<ve.size(); i++ ) {
156
157 ei = ve[i];
158
159 LOG("gcalchedisdiffxsec", pINFO) << "Energy: " << ei << " [GeV]";
160
161 TLorentzVector p4(0,0,ei,ei);
162 interaction->InitStatePtr()->SetProbeP4(p4);
163
164 for ( unsigned j=0; j<vy.size(); j++ ) {
165
166 double z = vy[j]; // z = (eo-em)/(ei-em)
167
168 by = (1-z)*(ei-10.)/ei; // y = 1 - eo/ei;
169
170 LOG("gcalchedisdiffxsec", pDEBUG) << " z: " << z << " y: " << by;
171
172 if ( by==1. ) by -= 1e-4;
173 else if ( by==0. ) {
174 by = 5./ei;
175 double by_prev = (1-vy[j-1])*(ei-10.)/ei;
176 LOG("gcalchedisdiffxsec", pDEBUG) << " by: " << by << " by_prev: " << by_prev << " " << vy[j-1];
177 if (by>by_prev) by = 1e-4;
178 }
179
180 if ( by<0 || by>1 ) continue;
181
182 interaction->KinePtr()->Sety(by);
183
184 if (fSaveAll) {
185 for ( unsigned k=0; k<vx.size(); k++ ) {
186 bx = vx[k];
187 interaction->KinePtr()->Setx(bx);
189 dxsec = xsec_alg->XSec(interaction, kPSxyfE) / units::cm2;
190 LOG("gcalchedisdiffxsec", pDEBUG) << "x: " << bx << " y: " << by << " -> d2sigmadxdy[E=" << ei << "GeV] = " << dxsec << " cm2";
191 tree->Fill();
192 }
193 }
194 else {
195 dxsec = 0.;
196 for ( unsigned k=0; k<vx.size(); k++ ) {
197 bx = vx[k];
198 interaction->KinePtr()->Setx(bx);
200 dxsec += xsec_alg->XSec(interaction, kPSxyfE) * bx;
201 }
202 dxsec *= widthx*TMath::Log(10) / units::cm2;
203 LOG("gcalchedisdiffxsec", pDEBUG) << " y: " << by << " -> d2sigmady[E=" << ei << "GeV] = " << dxsec << " cm2";
204 tree->Fill();
205 }
206
207 }
208
209 }
210
211 }
212
213 TFile * outfile = new TFile(fOutFileName.c_str(),"RECREATE");
214 tree->Write(tree->GetName());
215 outfile->Close();
216
217}
218
219
220//----------------------------------------------------------------------------------------
221//----------------------------------------------------------------------------------------
222//READ LIST
223//----------------------------------------------------------------------------------------
224//----------------------------------------------------------------------------------------
225vector<double> ReadListFromPath(string path) {
226
227 vector<double> list;
228
229 std::ifstream infile(path.c_str());
230
231 //check to see that the file was opened correctly:
232 if (!infile.is_open()) {
233 LOG("gcalchedisdiffxsec", pFATAL) << "There was a problem opening the input file!";
234 exit(1);//exit or do additional error checking
235 }
236
237 double val = 0.;
238 while (infile >> val) list.push_back(val);
239
240 return list;
241
242}
243
244
245//----------------------------------------------------------------------------------------
246//----------------------------------------------------------------------------------------
247//INPUT PARSER
248//----------------------------------------------------------------------------------------
249//----------------------------------------------------------------------------------------
250
251void DecodeCommandLine(int argc, char * argv[]) {
252
253 // Common run options. Set defaults and read.
255
256 // Parse run options for this app
257 CmdLnArgParser parser(argc,argv);
258
259 if( parser.OptionExists('p') ){
260 fNu = parser.ArgAsInt('p');
261 LOG("gcalchedisdiffxsec", pINFO) << "Probe = " << fNu;
262 }
263 else {
264 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified input neutrino type!";
265 PrintSyntax();
266 exit(1);
267 }
268
269 if( parser.OptionExists('t') ){
270 fTgt = parser.ArgAsInt('t');
271 LOG("gcalchedisdiffxsec", pINFO) << "Target = " << fTgt;
272 }
273 else {
274 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified input target type!";
275 PrintSyntax();
276 exit(1);
277 }
278
279 if( parser.OptionExists('x') ){
280 fPathXlist = parser.ArgAsString('x');
281 LOG("gcalchedisdiffxsec", pINFO) << "PathXlist = " << fPathXlist;
282 }
283 else {
284 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified input path to Xlist!";
285 PrintSyntax();
286 exit(1);
287 }
288
289 if( parser.OptionExists('y') ){
290 fPathYlist = parser.ArgAsString('y');
291 LOG("gcalchedisdiffxsec", pINFO) << "PathYlist = " << fPathYlist;
292 }
293 else {
294 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified input path to Ylist!";
295 PrintSyntax();
296 exit(1);
297 }
298
299 if( parser.OptionExists('e') ){
300 fPathElist = parser.ArgAsString('e');
301 LOG("gcalchedisdiffxsec", pINFO) << "PathElist = " << fPathElist;
302 }
303 else {
304 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified input path to Elist!";
305 PrintSyntax();
306 exit(1);
307 }
308
309 if( parser.OptionExists('s') ) {
310 fSaveAll = true;
311 LOG("gcalchedisdiffxsec", pINFO) << "SaveAll = " << fSaveAll;
312 }
313
314 if( parser.OptionExists('o') ){
315 fOutFileName = parser.ArgAsString('o');
316 LOG("gcalchedisdiffxsec", pINFO) << "OutFileName = " << fOutFileName;
317 }
318 else {
319 LOG("gcalchedisdiffxsec", pFATAL) << "Unspecified output file name!";
320 PrintSyntax();
321 exit(1);
322 }
323
324}
325
326//____________________________________________________________________________
327void PrintSyntax(void)
328{
329 LOG("gcalchedisdiffxsec", pNOTICE)
330 << "\n\n" << "Syntax:" << "\n"
331 << " gcalchedisdiffxsec -p nu -t tgt -o root_file\n"
332 << " -x pathxlist\n"
333 << " -y pathylist\n"
334 << " -e pathelist\n"
335 << " --tune genie_tune\n"
336 << " --event-generator-list list_name\n"
337 << " [--message-thresholds xml_file]\n";
338}
string infile
#define pNOTICE
Definition Messenger.h:61
#define pINFO
Definition Messenger.h:62
#define pFATAL
Definition Messenger.h:56
#define pDEBUG
Definition Messenger.h:63
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE,...
Definition Messenger.h:96
int main()
Command line argument parser.
string ArgAsString(char opt)
bool OptionExists(char opt)
was option set?
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
GENIE Event Generation Driver. A minimalist user interface object for generating neutrino interaction...
Definition GEVGDriver.h:54
void Configure(int nu_pdgc, int Z, int A)
const InteractionList * Interactions(void) const
void SetEventGeneratorList(string listname)
const EventGeneratorI * FindGenerator(const Interaction *interaction) const
Initial State information.
const Target & Tgt(void) const
void SetProbeP4(const TLorentzVector &P4)
A vector of Interaction objects.
Summary information for an interaction.
Definition Interaction.h:56
string AsString(void) const
const XclsTag & ExclTag(void) const
Definition Interaction.h:72
InitialState * InitStatePtr(void) const
Definition Interaction.h:74
const InitialState & InitState(void) const
Definition Interaction.h:69
Kinematics * KinePtr(void) const
Definition Interaction.h:76
void Setx(double x, bool selected=false)
void Sety(double y, bool selected=false)
string EventGeneratorList(void) const
Definition RunOpt.h:47
void ReadFromCommandLine(int argc, char **argv)
Definition RunOpt.cxx:99
void BuildTune()
build tune and inform XSecSplineList
Definition RunOpt.cxx:92
static RunOpt * Instance(void)
Definition RunOpt.cxx:54
int HitQrkPdg(void) const
Definition Target.cxx:242
bool HitSeaQrk(void) const
Definition Target.cxx:299
Cross Section Calculation Interface.
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
int FinalQuarkPdg(void) const
Definition XclsTag.h:72
string fOutFileName
string fPathYlist
void DecodeCommandLine(int argc, char *argv[])
string fPathElist
void PrintSyntax(void)
vector< double > ReadListFromPath(string path)
string fPathXlist
const double epsilon
bool fSaveAll
const double e
static constexpr double cm2
Definition Units.h:69
void MesgThresholds(string inpfile)
Definition AppInit.cxx:99
void UpdateWQ2FromXY(const Interaction *in)
THE MAIN GENIE PROJECT NAMESPACE
Definition AlgCmp.h:25