ROOT logo
 void testPDF(Double_t a = 208.)
{

    gSystem->Load("liblhapdf.so");
    StrucFunc_t pdf = kCTEQ6ll;
    

// Dynamically link some shared libs                    

    
//    
    char    parm[20][20];
    Double_t val[20];
//
//  Initialization
//
    printf("PDF initialized with: \n \n");
    strncpy(parm[0], "DEFAULT             ",  20);
    val[0]  = AliStructFuncType::PDFsetIndex(pdf);
    AliStructFuncType::PdfSet(parm, val);

//
//  Plots
//
// ================================================================
// Pdf 
//
    TCanvas *c1 = new TCanvas("c1","Gluon PDF",400,10,600,700);
    TF1* f_gl_2 = new TF1("f_gl_2", pdf_gl, -6., 0., 2);
    f_gl_2->SetParameter(0,2);
    f_gl_2->SetParameter(1,a);
    f_gl_2->SetMaximum(1000.);
    f_gl_2->SetLineColor(1);
    f_gl_2->Draw();

    TF1* f_gl_5 = new TF1("f_gl_5", pdf_gl, -6., 0., 2);
    f_gl_5->SetParameter(0,5);
    f_gl_5->SetParameter(1,a);
    f_gl_5->SetLineColor(2);
    f_gl_5->Draw("Same");

    TF1* f_gl_10 = new TF1("f_gl_10", pdf_gl, -6., 0., 2);
    f_gl_10->SetParameter(0,10);
    f_gl_10->SetParameter(1,a);
    f_gl_10->SetLineColor(3);
    f_gl_10->Draw("Same");

    TF1* f_gl_50 = new TF1("f_gl_50", pdf_gl, -6., 0., 2);
    f_gl_50->SetParameter(0,50);
    f_gl_50->SetParameter(1,a);
    f_gl_50->SetLineColor(4);
    f_gl_50->Draw("Same");

// ================================================================
// Gluon Modification
//
    TCanvas *c2 = new TCanvas("c2","Gluon PDF Nucl. Mod.",400,10,600,700);

    TF1* fmod_gl_2 = new TF1("fmod_gl_2", mod_gl, -6., 0., 2);
    fmod_gl_2->SetParameter(0,2.);
    fmod_gl_2->SetParameter(1,a);
    fmod_gl_2->SetMaximum(1.5);
    fmod_gl_2->SetLineColor(1);
    fmod_gl_2->Draw();
    fmod_gl_2->GetHistogram()->SetXTitle("log(x)");
    fmod_gl_2->GetHistogram()->SetYTitle("F^{g}_{A}/F^{g}_{p}");

    TF1* fmod_gl_5 = new TF1("fmod_gl_5", mod_gl, -6., 0., 2);
    fmod_gl_5->SetParameter(0,5.);
    fmod_gl_5->SetParameter(1,a);
    fmod_gl_5->SetLineColor(2);
    fmod_gl_5->Draw("Same");

    TF1* fmod_gl_10 = new TF1("fmod_gl_10", mod_gl, -6., 0., 2);
    fmod_gl_10->SetParameter(0,10.);
    fmod_gl_10->SetParameter(1,a);
    fmod_gl_10->SetLineColor(3);
    fmod_gl_10->Draw("Same");

    TF1* fmod_gl_50 = new TF1("fmod_gl_50", mod_gl, -6., 0., 2);
    fmod_gl_50->SetParameter(0,50.);
    fmod_gl_50->SetParameter(1,a);
    fmod_gl_50->SetLineColor(4);
    fmod_gl_50->Draw("Same");
    label();
    

//  =================================================================
//  Sea Quark Modification
//
    TCanvas *c3 = new TCanvas("c3","Sea Quark PDF Nucl. Mod.",400,10,600,700);
    TF1* fmod_sq_2 = new TF1("fmod_sq_2", mod_sq, -6., 0., 2);
    fmod_sq_2->SetParameter(0,2.);
    fmod_sq_2->SetParameter(1,a);
    fmod_sq_2->SetMaximum(1.5);
    fmod_sq_2->SetLineColor(1);
    fmod_sq_2->Draw();
    fmod_sq_2->GetHistogram()->SetXTitle("log(x)");
    fmod_sq_2->GetHistogram()->SetYTitle("F^{sea q}_{A}/F^{sea q}_{p}");

    TF1* fmod_sq_5 = new TF1("fmod_sq_5", mod_sq, -6., 0., 2);
    fmod_sq_5->SetParameter(0,5.);
    fmod_sq_5->SetParameter(1,a);
    fmod_sq_5->SetLineColor(2);
    fmod_sq_5->Draw("Same");

    TF1* fmod_sq_10 = new TF1("fmod_sq_10", mod_sq, -6., 0., 2);
    fmod_sq_10->SetParameter(0,10.);
    fmod_sq_10->SetParameter(1,a);
    fmod_sq_10->SetLineColor(3);
    fmod_sq_10->Draw("Same");

    TF1* fmod_sq_50 = new TF1("fmod_sq_50", mod_sq, -6., 0., 2);
    fmod_sq_50->SetParameter(0,50.);
    fmod_sq_50->SetParameter(1,a);
    fmod_sq_50->SetLineColor(4);
    fmod_sq_50->Draw("Same");
    label();
    
//  =================================================================
//  Valence Quark Modification
//
    TCanvas *c4 = new TCanvas("c4","Valence Quark PDF Nucl. Mod.",
			      400,10,600,700);
    TF1* fmod_vq_2 = new TF1("fmod_vq_2", mod_vq, -6., 0., 2);
    fmod_vq_2->SetParameter(0,2.);
    fmod_vq_2->SetParameter(1,a);
    fmod_vq_2->SetMaximum(1.5);
    fmod_vq_2->SetLineColor(1);
    fmod_vq_2->Draw();
    fmod_vq_2->GetHistogram()->SetXTitle("log(x)");
    fmod_vq_2->GetHistogram()->SetYTitle("F^{val q}_{A}/F^{val q}_{p}");

    TF1* fmod_vq_5 = new TF1("fmod_vq_5", mod_vq, -6., 0., 2);
    fmod_vq_5->SetParameter(0,5.);
    fmod_vq_5->SetParameter(1,a);
    fmod_vq_5->SetLineColor(2);
    fmod_vq_5->SetMaximum(1.5);
    fmod_vq_5->Draw("Same");

    TF1* fmod_vq_10 = new TF1("fmod_vq_10", mod_vq, -6., 0., 2);
    fmod_vq_10->SetParameter(0,10.);
    fmod_vq_10->SetParameter(1,a);
    fmod_vq_10->SetLineColor(3);
    fmod_vq_10->SetMaximum(1.5);
    fmod_vq_10->Draw("Same");

    TF1* fmod_vq_50 = new TF1("fmod_vq_50", mod_vq, -6., 0., 2);
    fmod_vq_50->SetParameter(0,50.);
    fmod_vq_50->SetParameter(1,a);
    fmod_vq_50->SetLineColor(4);
    fmod_vq_50->SetMaximum(1.5);
    fmod_vq_50->Draw("Same");
    label();
    
}
 

Double_t pdf_gl(Double_t* x, Double_t* par)
{
    Double_t xx = TMath::Power(10.,x[0]);
    Double_t y;
    Double_t upv, dnv, usea, dsea, str, chm, bot, top, gl;
    Double_t q = par[0];
    Double_t a = par[1];

    AliStructFuncType::StructA(xx, q, a, upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    y = gl;
    
    return y;
}

Double_t mod_gl(Double_t* x, Double_t* par)
{
    Double_t xx = TMath::Power(10.,x[0]);
    Double_t y;
    Double_t upv, dnv, usea, dsea, str, chm, bot, top, gl;
    Double_t q = par[0];
    Double_t a = par[1];

    AliStructFuncType::StructA(xx, q, a, upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    Double_t y1 = gl;

    AliStructFuncType::StructA(xx, q, 1., upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    return y1/gl;
}

Double_t mod_sq(Double_t* x, Double_t* par)
{
    Double_t xx = TMath::Power(10.,x[0]);
    Double_t y;
    Double_t upv, dnv, usea, dsea, str, chm, bot, top, gl;
    Double_t q = par[0];
    Double_t a = par[1];

    AliStructFuncType::StructA(xx, q, a, upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    Double_t y1 = usea+dsea;

    AliStructFuncType::StructA(xx, q, 1., upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    return y1/(usea+dsea);
}


Double_t mod_vq(Double_t* x, Double_t* par)
{
    Double_t xx = TMath::Power(10.,x[0]);
    Double_t y;
    Double_t upv, dnv, usea, dsea, str, chm, bot, top, gl;
    Double_t q = par[0];
    Double_t a = par[1];

    AliStructFuncType::StructA(xx, q, a, upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    Double_t y1 = upv+dnv;

    AliStructFuncType::StructA(xx, q, 1., upv, dnv, usea, 
				   dsea, str, chm, bot, top, gl);
    
    return y1/(upv+dnv);
}


void label()
{
    TLine *line = new TLine(-5.7,1.4,-5.2,1.4);
    line->SetLineColor(1);
    line->SetLineWidth(3);
    line->Draw();
    tex = new TLatex(-5.1, 1.38,"Q =   2 GeV");
    tex->SetTextSize(0.04);
    tex->SetLineWidth(2);
    tex->Draw();

    TLine *line = new TLine(-5.7,1.3,-5.2,1.3);
    line->SetLineColor(2);
    line->SetLineWidth(3);
    line->Draw();
    tex = new TLatex(-5.1, 1.28,"Q =   5 GeV");
    tex->SetTextSize(0.04);
    tex->SetLineWidth(2);
    tex->Draw();


    TLine *line = new TLine(-5.7,1.2,-5.2,1.2);
    line->SetLineColor(3);
    line->SetLineWidth(3);
    line->Draw();
    tex = new TLatex(-5.1, 1.18,"Q = 10 GeV");
    tex->SetTextSize(0.04);
    tex->SetLineWidth(2);
    tex->Draw();


    TLine *line = new TLine(-5.7,1.1,-5.2,1.1);
    line->SetLineColor(4);
    line->SetLineWidth(3);
    line->Draw();
    tex = new TLatex(-5.1, 1.08,"Q = 50 GeV");
    tex->SetTextSize(0.04);
    tex->SetLineWidth(2);
    tex->Draw();
}



 testPDF.C:1
 testPDF.C:2
 testPDF.C:3
 testPDF.C:4
 testPDF.C:5
 testPDF.C:6
 testPDF.C:7
 testPDF.C:8
 testPDF.C:9
 testPDF.C:10
 testPDF.C:11
 testPDF.C:12
 testPDF.C:13
 testPDF.C:14
 testPDF.C:15
 testPDF.C:16
 testPDF.C:17
 testPDF.C:18
 testPDF.C:19
 testPDF.C:20
 testPDF.C:21
 testPDF.C:22
 testPDF.C:23
 testPDF.C:24
 testPDF.C:25
 testPDF.C:26
 testPDF.C:27
 testPDF.C:28
 testPDF.C:29
 testPDF.C:30
 testPDF.C:31
 testPDF.C:32
 testPDF.C:33
 testPDF.C:34
 testPDF.C:35
 testPDF.C:36
 testPDF.C:37
 testPDF.C:38
 testPDF.C:39
 testPDF.C:40
 testPDF.C:41
 testPDF.C:42
 testPDF.C:43
 testPDF.C:44
 testPDF.C:45
 testPDF.C:46
 testPDF.C:47
 testPDF.C:48
 testPDF.C:49
 testPDF.C:50
 testPDF.C:51
 testPDF.C:52
 testPDF.C:53
 testPDF.C:54
 testPDF.C:55
 testPDF.C:56
 testPDF.C:57
 testPDF.C:58
 testPDF.C:59
 testPDF.C:60
 testPDF.C:61
 testPDF.C:62
 testPDF.C:63
 testPDF.C:64
 testPDF.C:65
 testPDF.C:66
 testPDF.C:67
 testPDF.C:68
 testPDF.C:69
 testPDF.C:70
 testPDF.C:71
 testPDF.C:72
 testPDF.C:73
 testPDF.C:74
 testPDF.C:75
 testPDF.C:76
 testPDF.C:77
 testPDF.C:78
 testPDF.C:79
 testPDF.C:80
 testPDF.C:81
 testPDF.C:82
 testPDF.C:83
 testPDF.C:84
 testPDF.C:85
 testPDF.C:86
 testPDF.C:87
 testPDF.C:88
 testPDF.C:89
 testPDF.C:90
 testPDF.C:91
 testPDF.C:92
 testPDF.C:93
 testPDF.C:94
 testPDF.C:95
 testPDF.C:96
 testPDF.C:97
 testPDF.C:98
 testPDF.C:99
 testPDF.C:100
 testPDF.C:101
 testPDF.C:102
 testPDF.C:103
 testPDF.C:104
 testPDF.C:105
 testPDF.C:106
 testPDF.C:107
 testPDF.C:108
 testPDF.C:109
 testPDF.C:110
 testPDF.C:111
 testPDF.C:112
 testPDF.C:113
 testPDF.C:114
 testPDF.C:115
 testPDF.C:116
 testPDF.C:117
 testPDF.C:118
 testPDF.C:119
 testPDF.C:120
 testPDF.C:121
 testPDF.C:122
 testPDF.C:123
 testPDF.C:124
 testPDF.C:125
 testPDF.C:126
 testPDF.C:127
 testPDF.C:128
 testPDF.C:129
 testPDF.C:130
 testPDF.C:131
 testPDF.C:132
 testPDF.C:133
 testPDF.C:134
 testPDF.C:135
 testPDF.C:136
 testPDF.C:137
 testPDF.C:138
 testPDF.C:139
 testPDF.C:140
 testPDF.C:141
 testPDF.C:142
 testPDF.C:143
 testPDF.C:144
 testPDF.C:145
 testPDF.C:146
 testPDF.C:147
 testPDF.C:148
 testPDF.C:149
 testPDF.C:150
 testPDF.C:151
 testPDF.C:152
 testPDF.C:153
 testPDF.C:154
 testPDF.C:155
 testPDF.C:156
 testPDF.C:157
 testPDF.C:158
 testPDF.C:159
 testPDF.C:160
 testPDF.C:161
 testPDF.C:162
 testPDF.C:163
 testPDF.C:164
 testPDF.C:165
 testPDF.C:166
 testPDF.C:167
 testPDF.C:168
 testPDF.C:169
 testPDF.C:170
 testPDF.C:171
 testPDF.C:172
 testPDF.C:173
 testPDF.C:174
 testPDF.C:175
 testPDF.C:176
 testPDF.C:177
 testPDF.C:178
 testPDF.C:179
 testPDF.C:180
 testPDF.C:181
 testPDF.C:182
 testPDF.C:183
 testPDF.C:184
 testPDF.C:185
 testPDF.C:186
 testPDF.C:187
 testPDF.C:188
 testPDF.C:189
 testPDF.C:190
 testPDF.C:191
 testPDF.C:192
 testPDF.C:193
 testPDF.C:194
 testPDF.C:195
 testPDF.C:196
 testPDF.C:197
 testPDF.C:198
 testPDF.C:199
 testPDF.C:200
 testPDF.C:201
 testPDF.C:202
 testPDF.C:203
 testPDF.C:204
 testPDF.C:205
 testPDF.C:206
 testPDF.C:207
 testPDF.C:208
 testPDF.C:209
 testPDF.C:210
 testPDF.C:211
 testPDF.C:212
 testPDF.C:213
 testPDF.C:214
 testPDF.C:215
 testPDF.C:216
 testPDF.C:217
 testPDF.C:218
 testPDF.C:219
 testPDF.C:220
 testPDF.C:221
 testPDF.C:222
 testPDF.C:223
 testPDF.C:224
 testPDF.C:225
 testPDF.C:226
 testPDF.C:227
 testPDF.C:228
 testPDF.C:229
 testPDF.C:230
 testPDF.C:231
 testPDF.C:232
 testPDF.C:233
 testPDF.C:234
 testPDF.C:235
 testPDF.C:236
 testPDF.C:237
 testPDF.C:238
 testPDF.C:239
 testPDF.C:240
 testPDF.C:241
 testPDF.C:242
 testPDF.C:243
 testPDF.C:244
 testPDF.C:245
 testPDF.C:246
 testPDF.C:247
 testPDF.C:248
 testPDF.C:249
 testPDF.C:250
 testPDF.C:251
 testPDF.C:252
 testPDF.C:253
 testPDF.C:254
 testPDF.C:255
 testPDF.C:256
 testPDF.C:257
 testPDF.C:258
 testPDF.C:259
 testPDF.C:260
 testPDF.C:261
 testPDF.C:262
 testPDF.C:263
 testPDF.C:264
 testPDF.C:265
 testPDF.C:266
 testPDF.C:267
 testPDF.C:268
 testPDF.C:269
 testPDF.C:270
 testPDF.C:271
 testPDF.C:272
 testPDF.C:273
 testPDF.C:274
 testPDF.C:275
 testPDF.C:276
 testPDF.C:277