Interface to CTEQ PDF Tables
Interface to CTEQ PDF Tables
v. 1.0.4
GENERAL DESCRIPTION
NEWS
New in the current version: This version is compatible with the new releases of the automake, autoconf and libtool packages. No changes in the code.
Current stable release: cteq-pdf-1.0.4.tar.gz
Previous stable release: cteq-pdf-1.0.3.tar.gz
Current development version:
Traditionally, the PDF functions are given in interpolation table form. Every PDF fits have their own table and a simple FORTRAN routine reads the table and does the interpolation. The drawback of this approach is that we cannot use more than one PDF table simultaneously. But sometimes it would be useful. For instance, we want to calculate the 3-jet cross sections at NLO for LHC. This is a hard and long calculation even with one PDF function; in some cases it requires days of CPU time. Of course we don’t want to repeat 40 times when we study the PDF uncertainties, we want to calculate the cross sections for the 40 PDFs in one run. To do this we need a new object oriented interface to the PDF tables.
I have developed new interfaces (C and F77) to the CTEQ PDF tables. These interfaces use the original CTEQ tables but the interface is implemented in object oriented way both in C and Fortran. This webpage is about this code.
DOWNLOAD AND INSTALL
Build from source code: In general this is the simplest thing to do. If you want to compile it from source you need C compiler. The base code doesn’t requires Fortran compile but make sure that you have a working Fortran compiler that is able to cooperate with your C compiler. I recommend to use GCC. It has C/C++ and F77/F95 compilers.
Current stable release: cteq-pdf-1.0.4.tar.gz
Previous stable release: cteq-pdf-1.0.3.tar.gz or Go to Archive
Current development version:
For Mac Users: If you have a Mac with MacOS X operating system you can install cteq-pdf using MacPorts system. You should do the following steps:
1. You have to install and activate the MacPorts system.
2. Edit the /opt/local/etc/ports/source.conf file and add the following line:
You can edit this file only with super user privilege.
3. In a shell run the following command:
sudo port -d sync
4. Now you can install the program by
sudo port install cteq-pdf
Now it is installed under the directory /opt/local/.
Now this installing procedure works with 1.8.2 version of MacPorts.
For Linux Users: A debian package would be nice....
DOCUMENTATIONS
Available PDF sets: All the current CTEQ PDF tables are hardwired in this release but you can use your any table file which is compatible with the original fortran interface. There are two way to refer to an internal table. You can use the Iset number (e.g.: Iset = 1 for CTEQ6M) or the name of the table (e.g.: CTEQ6M). I use the same Iset numbers names those are in the original FORTRAN interface. To get more detailed information about the available tables go to the CTEQ page or see Cteq6Pdf-2008.txt.
Test examples: There are some simple test example. Please check them to get some idea about basic concept.
Current test package: cteq-pdf-tests.zip
To run the test program you have to install the cteq-pdf library. Unpack the test package, edit the Makefile, make all and it should be ready to play with the simple programs (test1.c, test2.c test3.f test4.f test5.f test6.c).
Fortran Interface
At a certain level it is possible to do object oriented computing in Fortran. Of course the language doesn’t support it but we can code everything in C and make these functions callable in Fortran. Actually the whole Fortran interface is written in C. The main idea is to emulate pointers and this quasi pointers can refer to objects. Once we can refer to objects then we can define all sort of operations on them. We use the integer numbers as pointer for this purpose.
To create a PDF object and initialize with a PDF table one has to use either
function cteq_pdf_alloc_id(Iset)
integer Iset, cteq_pdf_alloc_id
or
function cteq_pdf_alloc_name(name)
integer cteq_pdf_alloc_name
character name(*)
The return value of these function is an integer and it refers to the allocated object. If it wasn’t able to create the object than it returns with value -1.
Example 1
.....
c----- A C/C++ programer never uses IMPLICIT declaration
implicit none
c----- Declare the external functions -----
integer cteq_pdf_alloc_id, cteq_pdf_alloc_name
c----- Declare the references to a PDF object -----
integer pdf(3)
c----- Create the objects ------
pdf(1) = cteq_pdf_alloc_id(100)
pdf(2) = cteq_pdf_alloc_name('CTEQ65.00')
pdf(3) = cteq_pdf_alloc_name('myTable.pds')
.....
In this example we have created three PDF objects. The first pdf(1) is initialized using the internal table files and the Iset number. Now pdf(1) refers to the CTEQ6M.00 table. The second object was initialized by using the name of the table and pdf(2) refers to the CTEQ65.00 table. You can also use your own table file which is not hard wired. If the table file is in the working directory (or you can provide the full path) then you can create an PDF object using the filename. This is in the third example. Now pdf(3) refers to your table. Note, the table file filename MUST have .pds or .tbl extension and MUST be compatible with the pds or tbl format, respectively. Of course, if you need more PDF object you can create more or you create then in a do loop.
Now, the objects are created and intialized. They are somewhere in the memory, who knows where, but we can refer to them. We want to use them. To calculate the PDF functions one has to use
function cteq_pdf_evolvepdf(pdf, iprtn, x, Q)
integer pdf, iprtn
double precision x, Q, cteq_pdf_alloc_id
Here the first argument is the reference to the PDF object. As we know it must be integer. The second argument is the parton label, iprtn = (6, 5, 4, 3, 2, 1, 0, -1, ..., -6) correspond to (t, b, c, s, d, u, g, u_bar, ..., t_bar). The last two arguments are the momentum fraction variable (x = 0,..1) and the scale parameter Q in GeV.
Example 2
.....
c----- Declare the external functions -----
double precision cteq_pdf_evolvepdf
c----- Declarations of local variables -----
integer iprtn
double precision x,q
c----- Plotting the gluon distributions -----
iprtn = 0
q = 1000.0d0
x = 0.1d0
10 write(*,*) x,(cteq_pdf_evolvepdf(pdf(i),iprtn,x,q), i=1,3)
x = x + 0.1d0
if(x < 1.0d0) goto 10
.....
In this example we use the previously allocated and initialized PDF tables simultaneously without re-initialization of the tables anytime when we switch from one table to another.
In the calculation we often need to compute the strong coupling. It is implemented in the library and can be called by using
function cteq_pdf_evolvas(pdf, Q)
integer pdf
double precision Q, cteq_pdf_alloc_id
Here the first argument is the reference to the PDF object. As we know it must be integer. The second argument is the scale parameter Q in GeV. This function returns with value of the strong coupling over 2π (return value is αs(Q)/2π).
The PDF object are allocated dynamically and we have free the memory if we don’t use them anymore. To do this use the following subroutine:
subroutine cteq_pdf_free(pdf)
integer pdf
Example 3
.....
c----- Free the pdf objects ------
do i=1,3
if(pdf(i).ne.-1) call cteq_pdf_free(pdf(i))
end do
.....
These are the most important functions. There are some other function to get some structural informations. Let us have a quick list of them.
function cteq_pdf_orderpd(pdf)
integer pdf, cteq_pdf_orderpdf
This returns the order of the PDF evolution. The return value is 1 for LO and 2 for NLO level evolution.
function cteq_pdf_orderas(pdf)
integer pdf, cteq_pdf_orderas
This returns the order of the strong coupling evolution. The return value is 1 for LO and 2 for NLO level evolution.
function cteq_pdf_nfmax(pdf)
integer pdf, cteq_pdf_orderas
This returns the maximum number of the active quark flavors.
function cteq_pdf_scale(pdf)
integer pdf
double precision cteq_pdf_scale
This returns the value of fitting scale in GeV.
function cteq_pdf_alfas(pdf)
integer pdf
double precision cteq_pdf_alfas
This returns the value of the strong coupling at the fitting scale.
function cteq_pdf_mass(pdf, iprtn)
integer pdf, iprtn
double precision cteq_pdf_mass
This returns the mass of the parton iprtn in GeV.
function cteq_pdf_threshold(pdf, nf)
integer pdf, nf
double precision cteq_pdf_threshold
This return the flavor thresholds in GeV.
Usage: If you use this library in your program then you have to link the cteq-pdf library when you compile your code. You should do something like
gfortran myProgram.f -o myProgram -L${HOME}/lib -lcteq-pdf
C Interface
With some propaganda purpose I started this documentation with the FORTRAN interface. The C documentation is coming. Meanwhile study the source code .......
C++ Interface
At the moment there is no C++ interface, but the C interface is fully object oriented and it shouldn’t be a problem to define a C++ class to wrap the C interface.