Skip to content

Framework for amplitude analysis in photoproduction processes

License

Notifications You must be signed in to change notification settings

dwinney/jpacPhoto

 
 

Folders and files

NameName
Last commit message
Last commit date

Latest commit

 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 

Repository files navigation

jpacPhoto

Framework for amplitude analysis involving single meson production via quasi-elastic scattering on a nucleon target. Focus on expandability and easy interfacing with other libraries / analysis code.

Compilation of the base library requires only ROOT (tested with versions 6.17 and 6.24) with MathMore and Boost C++ (version $\geq$ 1.68) libraries.

INSTALLATION

To install, clone normally and use:

cd jpacPhoto
mkdir build && cd build
cmake ..
cmake --build . --target install

This will create the core library /lib/libJPACPHOTO.so with the linkable library as well as ROOT dictionary (.pcm) files.

Additionally a scripting executable will be installed into /bin/jpacPhoto which short-cuts loading the libraries into an interactive ROOT session and running a .cpp file as an interpreted script. This executable requires the environment variable JPACPHOTO to be set to the top-level directory in order to find auxilary files. This can be done as such:

export JPACPHOTO=/path/to/jpacPhoto # for bash
setenv JPACPHOTO /path/to/jpacPhoto # for csh

USAGE

The primary use case is to reproduce results from recent [1-5] and older analyses [6-8] from the JPAC collaboration in a general and unified framework for amplitude analysis in photoproduction. The compiled library contains the framework of amplitudes and observables as abstract classes with specific amplitude models needed to be imported at run-time as header-only libraries.

The compiled jpacPhoto executable pipes an analysis script, relevent amplitude files, and the compiled library into ROOT's cling interpeter to run. This set up mimics a Python-like environment without requiring recompilation of the whole library when changes are made to amplitude files. Amplitudes and scripts relevant for JPAC papers are located in /physics and /scripts respectively. To run a script located in the bin directory simply run

jpacPhoto my_script.cpp

or add the bin directory to $PATH to call jpacPhoto from any directory.

The linkable library can be added to an existing CMake project using find_library() and linked as normal:

find_library(JPACPHOTO NAMES JPACPHOTO libJPACPHOTO
                       HINTS "$ENV{JPACPHOTO}/lib")
target_link_libraries( myTarget JPACPHOTO)

AMPLITUDES

The main object of interest in the core library is the abstract amplitude and derived implementations defined by the user for specific physics models. All amplitudes available so far may be found in /physics as well as a template file to guide adding new classes. Amplitudes are calculated on a per-helicity basis which allows one to compute an array of observables (units of GeV and nb assumed where appropriate):

Observable Callable amplitude function
Unpolarized probability distribution $\sum_{{\lambda}} |A_{{\lambda}}|^2$ probability_distribution(double s, double t)
Differential cross section $d\sigma/dt$ differential_xsection(double s, double t)
Polarized differential cross section $d\sigma_{\perp,\parallel}/dt$ polarized_dxsection(double pm, double s, double t)
Integrated total cross section $\sigma$ integrated_xsection(double s)
Double polarization asymmetries $A_{LL},K_{LL}$ A_LL(double s, double t)
K_LL(double s, double t)
Meson spin density matrix elements $\rho^{\alpha}_{\lambda,\lambda^\prime}$ mSDME_H(int a, int lam, int lamp, double s, double t)
mSDME_GJ(int a, int lam, int lamp, double s, double t)
Baryon spin density matrix elements $\rho^{\alpha}_{\frac{\lambda}{2},\frac{\lambda^\prime}{2}}$ bSDME_H(int a, int lam, int lamp, double s, double t)
bSDME_GJ(int a, int lam, int lamp, double s, double t)
Beam asymmetry $\Sigma_{4\pi}$ beam_asymmetry_4pi(double s, double t)

All kinematics are passed around by the kinematics class which allows arbitrary masses for all particles and arbitrary quantum numbers for the produced final state meson & baryon. Multiple amplitudes may describe the same process and share the same kinematics instance.

The basic usage is:

// Set up kinematics
kinematics myKin = new_kinematics(M_MESON, M_BARYON);
myKin->set_meson_JP(1, -1);  // J = 1  , P = -1
myKin->set_baryon_JP(1, 1);  // J = 1/2, P =  1

// Set up amplitude
amplitude myAmp1 = new_amplitude<my_implementation>(myKin, /* additional parameters */);
myAmp1->set_parameters{ {/* couplings etc */} };

amplitude myAmp2 = new_amplitude<my_other_implementation>(myKin, /* additional parameters */);
myAmp2->set_parameters{ {/* couplings etc */} };

// Evaluate observables
myAmp1->integrated_xsection(s);
myAmp2->mSDME_GJ(0, 1, -1, s, t);

Interfering sums of amplitues and/or partial wave projections can be constucted out of any combination of compatible amplitude instances. These are all treated as amplitudes themselves and have access to all observables with the same syntax:

// Two different amplitudes but share same kinematics
amplitude amp1, amp2;

// Sum together
amplitude sum = amp1 + amp2;
sum->set_parameters{ /* couplings for both amp1 and amp2 */};

// Take the J = 1 partial wave
amplitude pwave = project(1, sum);

// Access observables
amp1->integrated_xsection(s);  // Individual term
sum->integrated_xsection(s);   // Interfering sum
pwave->integrated_xsection(s); // Only P-wave contribution of sum

SEMI-INCLUSIVE DISTRIBUTIONS

Methods for semi-inclusive processes can be added via the semi_inclusive class. These are used for example in [3-4] to investigate inclusive XYZ production. Because semi-inclusive models are implemented at the cross section level and lose explicit helicity dependence, only unpolarized observables are available (Units of GeV and nb assumed where appropriate):

Observable Callable semi_inclusive function
Lorentz-invariant cross section $E_{\mathcal{Q}} \frac{d^3\sigma}{d^3\mathbf{q}}$ invariant_xsection(double s, double t, double M2)
Double-differential cross section $\frac{d^2\sigma}{dt dM^2}$,
$\frac{d^2\sigma}{dt dx}$,
$\frac{d^2\sigma}{dx dy^2}$,
dsigma_dtdM2(double s, double t, double M2)
dsigma_dtdx(double s, double t, double x)
dsigma_dxdy2(double s, double x, double y2)
Single-differential cross section $\frac{d\sigma}{dt}$,
$\frac{d\sigma}{dM^2}$,
$\frac{d\sigma}{dx}$,
$\dots$
dsigma_dt(double s, double t)
dsigma_dM2(double s, double M2)
dsigma_dx(double s, double x)
$\dots$
Fully integrated cross section $\sigma$ integrated_xsection(double s)

Much of the syntax is the same as with exclusive amplitudes although models are implemented at the level of amplitude-squared.

// Use same interface with amplitudes
kinematics kZ = new_kinematics(M_Z);

// Semi-inclusive Z production
semi_inclusive Z_inc = new_semi_inclusive<inclusive_model>(kZ, /* other parameter */);
Z_inc->set_parameters({ /* couplings, etc */ });

// Some processes we want to add the exclusive amplitude
amplitude Z_exc = new_amplitude<exclusive_model>(kZ, /* etc */);
Z_inc += Z_exc; // Add to the purely inclusive distribution

// Get observables
Z_inc->integrated_xsection(s);

ANALYSIS TOOLS

Tools to fit amplitudes to experimental data are available through the fitter and plotter classes. Data may be imported using the data_set class as interface. Arbitrarily many data sets may be imported into a fitter where one must specify the minimazition function per data type. An end-to-end example used in [4] may be found in the appropriate scripts directory.

A schematic analysis may look like:

struct my_analysis
{
    // Strings to identify data types in fitter's terminal messages
    static std::string data_type(int i)
    {
        if (i == 0) "Allowed type";
        else        "Not allowed type";
    };

    // Minimization funcion
    static double fcn(const std::vector<data_set> &data, amplitude amp)
    {
        double chi2 = 0;
        for (auto datum : data) if (data._type == 0) chi2 += /* calculate chi2*/;
        return chi2;
    };
};

int main()
{
  // Set up amplitude
  amplitude my_amplitude;

  // Do fit
  fitter fitter<my_analysis>(my_amplitude);
  fitter.add_data(my_data);
  fitter.set_parameter_limits("parameter 1", {0., 1});
  fitter.fix_parameter("fixed parameter", 0.5);
  fitter.sync_parameter("parameter 2", "parameter 3"); // enforce 2 == 3 with 3 free
  fitter.do_fit();

  // Plot results
  plotter plotter;
  plot p = plotter.new_plot();
  p.add_data(my_data);
  p.add_curve({min, max}, [&](double x){ return my_amplitude->observable(x); }, "label");
  p.save("outfile.pdf");
};

AmpTools

In addition to the built-in analysis tools, the aim of this project is to be usable with existing tools. A relatively simple interface with AmpTools is provided in /AmpTools as well as scripts in the analogous /scripts directory. More complicated processes involving unstable particles can be incorporated by wrapping the jpacPhoto helicity amplitudes into AmpTools classes which decay lineshapes.

To run these scripts with the jpacPhoto executable requires having AmpTools to be built and the AMPTOOLS environment variable set to the top level install directory. Since the cling interpreter can only load dynamic libraries (which is not the default installation mode of AmpTools), the optional cmake flag -DDYNAMIC_AMPTOOLS=TRUE is provided to repackage an existic static library to one shared version.

REFERENCES

About

Framework for amplitude analysis in photoproduction processes

Resources

License

Stars

Watchers

Forks

Releases

No releases published

Packages

No packages published

Languages

  • C++ 98.9%
  • Other 1.1%