Page MenuHomeHEPForge

tableGeneratorMain.cpp
No OneTemporary

tableGeneratorMain.cpp

//==============================================================================
// tableGeneratorMain.cpp
//
// Copyright (C) 2010-2013 Tobias Toll and Thomas Ullrich
//
// This file is part of Sartre version: 1.1
//
// This program is free software: you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation.
// This program is distributed in the hope that it will be useful,
// but without any warranty; without even the implied warranty of
// merchantability or fitness for a particular purpose. See the
// GNU General Public License for more details.
// You should have received a copy of the GNU General Public License
// along with this program. If not, see <http://www.gnu.org/licenses/>.
//
// Author: Thomas Ullrich
// Last update:
// $Date: 2013-05-29 21:27:49 +0100 (Wed, 29 May 2013) $
// $Author: thomas.ullrich@bnl.gov $
//==============================================================================
//
// Main program to create amplitude lookup tables.
// [Developer only]
//==============================================================================
#include <sstream>
#include <cstdlib>
#include <iostream>
#include "Amplitudes.h"
#include "TH1D.h"
#include "TFile.h"
#include "Constants.h"
#include "Table.h"
#include "TableGeneratorSettings.h"
#include "Enumerations.h"
#define PR(x) cout << #x << " = " << (x) << endl;
using namespace std;
int main(int argc, char *argv[]){
//To explicitly delete all new histograms myself:
TH1::AddDirectory(false);
TableGeneratorSettings* settings = TableGeneratorSettings::instance();
//
// Check arguments
//
char* runcard;
if (argc != 4) {
cout<<"Usage: tableGeneratorMain runcard startBin endBin"<<endl;
return 1;
}
else{
runcard = argv[1];
settings->setStartBin(atoi(argv[2]));
settings->setEndBin(atoi(argv[3]));
}
settings->readSettingsFromFile(runcard);
settings->consolidateSettings();
int nBinQ2 = settings->Q2bins();
int nBinW2 = settings->W2bins();
int nBinT= settings->tbins();
double Q2min= settings->Q2min();
double Q2max= settings->Q2max();
double Wmin= settings->Wmin();
double Wmax= settings->Wmax();
double W2min=Wmin*Wmin;
double W2max=Wmax*Wmax;
double tmin= settings->tmin();
double tmax= settings->tmax();
unsigned int massA=settings->A();
int vmPDG = settings->vectorMesonId();
DipoleModelType model = settings->dipoleModel();
int startingBin=settings->startBin();
int endingBin=settings->endBin();
int modes=settings->modesToCalculate();
Table tableT;
Table tableL;
Table tableT2;
Table tableL2;
bool logQ2=true, logW2=false, logT=false, logC=true;
//set filenames for the tables:
string rootfile=settings->rootfile();
ostringstream filenameT, filenameL, filenameT2, filenameL2;
filenameT.str("");
filenameT<<rootfile<<"_bin"<<startingBin<<"_T.root";
filenameL.str("");
filenameL<<rootfile<<"_bin"<<startingBin<<"_L.root";
filenameT2.str("");
filenameT2<<rootfile<<"_bin"<<startingBin<<"_T2.root";
filenameL2.str("");
filenameL2<<rootfile<<"_bin"<<startingBin<<"_L2.root";
(void) tableT.create(nBinQ2, Q2min, Q2max,
nBinW2, W2min, W2max,
nBinT, tmin, tmax,
logQ2, logW2, logT, logC, // all bools
mean_A, transverse,
massA, vmPDG, model,
filenameT.str().c_str());
(void) tableL.create(nBinQ2, Q2min, Q2max,
nBinW2, W2min, W2max,
nBinT, tmin, tmax,
logQ2, logW2, logT, logC, // all bools
mean_A, longitudinal,
massA, vmPDG, model,
filenameL.str().c_str());
(void) tableT2.create(nBinQ2, Q2min, Q2max,
nBinW2, W2min, W2max,
nBinT, tmin, tmax,
logQ2, logW2, logT, logC, // all bools
mean_A2, transverse,
massA, vmPDG, model,
filenameT2.str().c_str());
(void) tableL2.create(nBinQ2, Q2min, Q2max,
nBinW2, W2min, W2max,
nBinT, tmin, tmax,
logQ2, logW2, logT, logC, // all bools
mean_A2, longitudinal,
massA, vmPDG, model,
filenameL2.str().c_str());
//Create and initialize the amplitudes calculator:
Amplitudes amps;
//Generate the the nucleon configurations:
amps.generateConfigurations();
// Print out settings:
cout<<"Tables will be generated for:"<<endl;
cout<<"Mode to calculate: "<<modes<<endl;
cout<<"Vector Meson Id: "<<vmPDG<<endl;
cout<<"Bins: "<<startingBin<<"-"<<endingBin<<endl;
cout<<"Q2 range: ["<<Q2min<<", "<<Q2max<<"], "<<nBinQ2<<" bins."<<endl;
cout<<"W2 range: ["<<W2min<<", "<<W2max<<"], "<<nBinW2<<" bins."<<endl;
cout<<" t range: ["<<tmin<<", "<<tmax<<"], "<<nBinT<<" bins."<<endl;
cout<<"Dipole Model: "<<model<<endl;
cout<<endl;
cout<<"Filling tables..."<<endl;
time_t tableFillStart = time(0);
for(int i=startingBin; i<endingBin; i++) {
double Q2, W2, t;
tableT.binCenter(i, Q2, W2, t);
//calculate contents and fill tables:
amps.calculate(t, Q2, W2);
double aT=amps.amplitudeT();
double aL=amps.amplitudeL();
if(massA>1){
tableT.fill(i, aT);
tableL.fill(i, aL);
}
if(modes!=1){
double aT2=amps.amplitudeT2();
double aL2=amps.amplitudeL2();
tableT2.fill(i, aT2);
tableL2.fill(i, aL2);
}
}
cout<<"CPU Time/Entry: " << double(time(0)-tableFillStart)/endingBin<<" s/entry."<<endl;
cout<<"Total time: "<<double(time(0)-tableFillStart)/60./60.<<" h"<<endl;
if(massA>1){
tableT.write();
tableL.write();
}
if(modes!=1){
tableT2.write();
tableL2.write();
}
return 0;
}

File Metadata

Mime Type
text/x-c
Expires
Sat, Dec 21, 4:49 PM (18 h, 45 m)
Storage Engine
blob
Storage Format
Raw Data
Storage Handle
4014260
Default Alt Text
tableGeneratorMain.cpp (6 KB)

Event Timeline