BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
#ifndef BRAT_BrZdcPedCalModule
#include "BrZdcPedCalModule.h"
#endif

#ifndef BRAT_BrTableNames
#include "BrTableNames.h"
#endif

#ifndef BRAT_BrZdcDig
#include "BrZdcDig.h"
#endif

#ifndef ROOT_TF1
#include "TF1.h"
#endif

ClassImp (BrZdcPedCalModule);

 BrZdcPedCalModule::BrZdcPedCalModule ()
{
    SetState (kSetup);
}



 BrZdcPedCalModule::BrZdcPedCalModule (const Char_t *name, const Char_t *title):BrZdcCalModule (name, title)
{
    SetState (kSetup);
}



 BrZdcPedCalModule::~BrZdcPedCalModule ()
{

}



 void BrZdcPedCalModule::DefineHistograms ()
{
    if (fCommitAscii || fLoadAscii) return;
    if (GetState ()!= kInit) return;

    TDirectory *saveDir= gDirectory;
    TDirectory* histDir= gDirectory->mkdir (Form (GetName ()));

    histDir->cd ();

    fPedHi [0]= new TH1F ("ped1", "", 100, 0., 100.);
    fPedHi [1]= new TH1F ("ped2", "", 100, 0., 100.);
    fPedHi [2]= new TH1F ("ped3", "", 100, 0., 100.);
    fPedHi [3]= new TH1F ("pedSum", "", 100, 0., 100.);

    fPedLo [0]= new TH1F ("ped1_Lo", "", 100, 0., 100.);
    fPedLo [1]= new TH1F ("ped2_Lo", "", 100, 0., 100.);
    fPedLo [2]= new TH1F ("ped3_Lo", "", 100, 0., 100.);
    fPedLo [3]= new TH1F ("pedSum_Lo", "", 100, 0., 100.);

    gDirectory= saveDir;
}



 void BrZdcPedCalModule::Init ()
{
    SetState (kInit);

    BrZdcCalModule::Init ();
    BrCalibration::EAccessMode mode= BrCalibration::kRead;

    if (fLoadAscii) mode= BrCalibration::kTransitional;
    else if (fSaveAscii || fCommitAscii) mode= BrCalibration::kWrite;

    fCalibration->Use ("PedHi", mode, fNoModules);
    fCalibration->Use ("PedWidthHi", mode, fNoModules);
    fCalibration->Use ("PedLo", mode, fNoModules);
    fCalibration->Use ("PedWidthLo", mode, fNoModules);
}



 void BrZdcPedCalModule::Begin ()
{
    SetState (kBegin);

    if (fLoadAscii) {ReadAscii (); return;}
    if (fCommitAscii) return;

    for (Int_t i= 0; i < fNoModules; i++)
    {
        fCalibration->SetPedHi (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetPedWidthHi (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetPedLo (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetPedWidthLo (i+1, BrZdcCalibration::kCalException);
    }
}



 void BrZdcPedCalModule::Event (BrEventNode *inNode, BrEventNode *outNode)
{
    SetState (kEvent);

    if (fCommitAscii || fLoadAscii) return;

    BrZdcDig *digzdc_p;

    if ((digzdc_p= (BrZdcDig*) inNode->GetObject (kDigZDC)))
    {
        if (!strcmp (GetName (), "ZDCLeft"))
        {
            fPedHi [0]->Fill (digzdc_p->GetLeftAdc (0));
            fPedHi [1]->Fill (digzdc_p->GetLeftAdc (1));
            fPedHi [2] ->Fill (digzdc_p->GetLeftAdc (2));
            fPedHi [3]->Fill (digzdc_p->GetLeftAdcSum ());
        }
        if (!strcmp (GetName (), "ZDCRight"))
        {
            fPedHi [0]->Fill (digzdc_p->GetRightAdc (0));
            fPedHi [1]->Fill (digzdc_p->GetRightAdc (1));
            fPedHi [2]->Fill (digzdc_p->GetRightAdc (2));
            fPedHi [3]->Fill (digzdc_p->GetRightAdcSum ());
        }
    }

    if ((digzdc_p= (BrZdcDig*) inNode->GetObject (kDigZDCLo)))
    {
        if (!strcmp (GetName (), "ZDCLeft"))
        {
            fPedLo [0]->Fill (digzdc_p->GetLeftAdc (0));
            fPedLo [1]->Fill (digzdc_p->GetLeftAdc (1));
            fPedLo [2] ->Fill (digzdc_p->GetLeftAdc (2));
            fPedLo [3]->Fill (digzdc_p->GetLeftAdcSum ());
        }
        if (!strcmp (GetName (), "ZDCRight"))
        {
            fPedLo [0]->Fill (digzdc_p->GetRightAdc (0));
            fPedLo [1]->Fill (digzdc_p->GetRightAdc (1));
            fPedLo [2]->Fill (digzdc_p->GetRightAdc (2));
            fPedLo [3]->Fill (digzdc_p->GetRightAdcSum ());
        }
    }
}



 void BrZdcPedCalModule::End ()
{
    SetState (kEnd);
}



 void BrZdcPedCalModule::Finish ()
{
    SetState (kFinish);

    if (fLoadAscii) return;
    if (fCommitAscii) {ReadAscii (); return;}

    Int_t i;

    Float_t min1, min2, min3, minSum;
    Float_t max1, max2, max3, maxSum;
    Float_t meanHi [4], sigmaHi [4], meanLo [4], sigmaLo [4];

    TF1 *pedFuncHi [4], *pedFuncLo [4];

    //----- high-gain ADC pedestals

    min1= (fPedHi [0]->GetMaximumBin ()- 5);
    min2= (fPedHi [1]->GetMaximumBin ()- 5);
    min3= (fPedHi [2]->GetMaximumBin ()- 5);
    minSum= (fPedHi [3]->GetMaximumBin ()- 5);

    max1= (fPedHi [0]->GetMaximumBin ()+ 5);
    max2= (fPedHi [1]->GetMaximumBin ()+ 5);
    max3= (fPedHi [2]->GetMaximumBin ()+ 5);
    maxSum= (fPedHi [3]->GetMaximumBin ()+ 5);

    pedFuncHi [0]= new TF1 ("func1", "gaus", min1, max1);
    pedFuncHi [1]= new TF1 ("func2", "gaus", min2, max2);
    pedFuncHi [2]= new TF1 ("func3", "gaus", min3, max3);
    pedFuncHi [3]= new TF1 ("funcSum", "gaus", minSum, maxSum);

    fPedHi [0]->Fit ("func1", "QR");
    fPedHi [1]->Fit ("func2", "QR");
    fPedHi [2]->Fit ("func3", "QR");
    fPedHi [3]->Fit ("funcSum", "QR");

    for (i= 0; i < fNoModules; i++)
    {
        meanHi [i]= pedFuncHi [i]->GetParameter (1);
        sigmaHi [i]= pedFuncHi [i]->GetParameter (2);
    }

    //----- low-gain ADC pedestals

    min1= (fPedLo [0]->GetMaximumBin ()- 5);
    min2= (fPedLo [1]->GetMaximumBin ()- 5);
    min3= (fPedLo [2]->GetMaximumBin ()- 5);
    minSum= (fPedLo [3]->GetMaximumBin ()- 5);

    max1= (fPedLo [0]->GetMaximumBin ()+ 5);
    max2= (fPedLo [1]->GetMaximumBin ()+ 5);
    max3= (fPedLo [2]->GetMaximumBin ()+ 5);
    maxSum= (fPedLo [3]->GetMaximumBin ()+ 5);

    pedFuncLo [0]= new TF1 ("func1_Lo", "gaus", min1, max1);
    pedFuncLo [1]= new TF1 ("func2_Lo", "gaus", min2, max2);
    pedFuncLo [2]= new TF1 ("func3_Lo", "gaus", min3, max3);
    pedFuncLo [3]= new TF1 ("funcSum_Lo", "gaus", minSum, maxSum);

    fPedLo [0]->Fit ("func1_Lo", "QR");
    fPedLo [1]->Fit ("func2_Lo", "QR");
    fPedLo [2]->Fit ("func3_Lo", "QR");
    fPedLo [3]->Fit ("funcSum_Lo", "QR");


    for (i= 0; i < fNoModules; i++)
    {
        meanLo [i]= pedFuncLo [i]->GetParameter (1);
        sigmaLo [i]= pedFuncLo [i]->GetParameter (2);
    }

    //----- set parameters for eventual commission into DB

    for (i= 0; i < fNoModules; i++)
    {
        fCalibration->SetPedHi (i+1, meanHi [i]);
        fCalibration->SetPedWidthHi (i+1, sigmaHi [i]);
        fCalibration->SetPedLo (i+1, meanLo [i]);
        fCalibration->SetPedWidthLo (i+1, sigmaLo [i]);
    }

    fCalibration->SetComment ("PedHi", "Generated by  BrZdcPedCalModule: fit with a gaussian function");
    fCalibration->SetComment ("PedWidthHi", "Generated by  BrZdcPedCalModule: fit with a gaussian function");
    fCalibration->SetComment ("PedLo", "Generated by  BrZdcPedCalModule: fit with a gaussian function");
    fCalibration->SetComment ("PedWidthLo", "Generated by  BrZdcPedCalModule: fit with a gaussian function");

    if (fSaveAscii) SaveAscii ();
}



 void BrZdcPedCalModule::SaveAscii ()
{
    //----- save pedestals to an ASCII file

    Int_t i;
    Float_t par1, par2, par3, par4;

    BrRunInfoManager *runMan= BrRunInfoManager::Instance ();
    const BrRunInfo *run= runMan->GetCurrentRun ();

    if (run->GetRunNo ()== -1) { Abort ("SaveAscii", "RunInfo has run number = -1"); return; }

    BrZdcCalModule::SaveAscii ();
    ofstream file (fCalibFile.Data (), ios::out);

    file.setf (ios::left, ios::adjustfield);
    file.setf (ios::showpoint);
    file.precision (3);
    file.setf (ios::fixed, ios::floatfield);

    file << "**************************************************" << endl;
    file << "*" << endl;
    file << "*" << "	" << GetName () << " pedestals: events from run #" << run->GetRunNo () << endl;
    file << "*" << endl;
    file << "**************************************************" << endl;
    file << "*" << endl;
    file << "*" << "	";

    file.width (8);
    file << "Mod.#";

    file.width (15);
    file << "MeanPedHi";

    file.width (15);
    file << "SigmaPedHi";

    file.width (15);
    file << "MeanPedLo";

    file.width (15);
    file << "SigmaPedLo" << endl;

    file << "*" << endl;

    for (i= 0; i < fNoModules; i++)
    {
        par1= fCalibration->GetPedHi (i+1);
        par2= fCalibration->GetPedWidthHi (i+1);
        par3= fCalibration->GetPedLo (i+1);
        par4= fCalibration->GetPedWidthLo (i+1);

        file << "	";
        file.width (8);
        file << (i+1);

        file.width (15);
        file << par1;

        file.width (15);
        file << par2;

        file.width (15);
        file << par3;

        file.width (15);
        file << par4 << endl;
    }
}



 void BrZdcPedCalModule::ReadAscii ()
{
    //----- save pedestals to an ASCII file

    BrZdcCalModule::SaveAscii ();
    ifstream file (fCalibFile.Data (), ios::in);

    if (!file) { Abort ("ReadAscii", "File %s was not found", fCalibFile.Data ()); return; }

    Int_t i= 0, module;
    Float_t par1, par2, par3, par4;
    Char_t c;

    do
    {
        file.get (c);
        if (c== '*') file.ignore (256, 'n');
    }
    while (c== '*');

    for (Int_t i= 0; i < fNoModules; i++)
    {
        file >> module >> par1 >> par2 >> par3 >> par4;

        fCalibration->SetPedHi (i+1, par1);
        fCalibration->SetPedWidthHi (i+1, par2);
        fCalibration->SetPedLo (i+1, par3);
        fCalibration->SetPedWidthLo (i+1, par4);
    }
}

This page automatically generated by script docBrat by Christian Holm

Copyright ; 2002 BRAHMS Collaboration <brahmlib@rcf.rhic.bnl.gov>
Last Update on by

Validate HTML
Validate CSS