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

#ifndef BRAT_BrEvent
#include "BrEvent.h"
#endif

#ifndef BRAT_BrEventHeader
#include "BrEventHeader.h"
#endif

#ifndef BRAT_BrBbVertex
#include "BrBbVertex.h"
#endif

#ifndef BRAT_BrDataObject
#include "BrDataObject.h"
#endif

#ifndef BRAT_BrZdcRdo
#include "BrZdcRdo.h"
#endif

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

#ifndef ROOT_TH1
#include "TH1.h"
#endif

#ifndef ROOT_TNtuple
#include "TNtuple.h"
#endif


ClassImp (BrZdcSlewCalModule);

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



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



 BrZdcSlewCalModule::~BrZdcSlewCalModule ()
{

}



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

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

    histDir->cd ();

    fSlewProf [0]= new TProfile ("slew1", "module1", 200, 0., 0.2);
    fSlewProf [1]= new TProfile ("slew2", "module2", 200, 0., 0.2);
    fSlewProf [2]= new TProfile ("slew3", "module3", 200, 0., 0.2);
    fSlewProf [3]= new TProfile ("slewSum", "sum", 200, 0., 0.2);

    fTime [0]= new TH1F ("time1", "module1", 100, 0., 100.0);
    fTime [1]= new TH1F ("time2", "module2", 100, 0., 100.0);
    fTime [2]= new TH1F ("time3", "module3", 100, 0., 100.0);

    if(fUseNtuple)
      fNtuple= new TNtuple ("BBZDC", "BB-ZDC ntuple",
			    "bbLeftTime:bbRightTime:bbZ0:bbTime0:bbMethod:"
			    "zdcLeftTime1:zdcLeftTime2:zdcLeftTime3:zdcRightTime1:zdcRightTime2:zdcRightTime3:"
			    "leftAdc1:leftAdc2:leftAdc3:rightAdc1:rightAdc2:rightAdc3");
    
    gDirectory= saveDir;
}



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

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

    fCalibration->Use ("TdcGain", mode, fNoModules); // use TdcGain for slew-correction

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

    fCalibration->Use ("Slewpar1", mode, fNoModules);
    fCalibration->Use ("Slewpar2", mode, fNoModules);
    fCalibration->Use ("Slewpar3", mode, fNoModules);
    fCalibration->Use ("Slewpar4", mode, fNoModules);
    fCalibration->Use ("Slewpar5", mode, fNoModules);
}



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

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

    for (Int_t i= 0; i < fNoModules; i++)
    {
        fCalibration->SetSlewpar1 (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetSlewpar2 (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetSlewpar3 (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetSlewpar4 (i+1, BrZdcCalibration::kCalException);
        fCalibration->SetSlewpar5 (i+1, BrZdcCalibration::kCalException);
    }
}



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

    if (fCommitAscii || fLoadAscii) return;
    if (inNode->IsA ()!= BrEvent::Class ())
    {
        Failure ("Event", "input node not an event");
        return;
    }

    Int_t i, trig, trig1, bbMethod;
    Float_t bbTime= 0.0, adc [4], tdc [4], tdcGain [4];

    for (i= 0; i < fNoModules; i++) { adc [i]= 0.0; tdc [i]= 0.0; }

    BrEvent *e= (BrEvent*) inNode;
    BrEventHeader *header= e->GetEventHeader ();

    trig= header->TriggerWord (1);
    trig1= (0x0008 & trig); // ZDC trigger before downscaling

    if (trig1)
    {
        BrDataObject *bb= 0;
        bb= inNode->GetObject ("BB Vertex");
        if (!bb) return;

        bbMethod= ((BrBbVertex*) bb)->GetMethod ();
        if (bbMethod!= 2) return;

        BrZdcRdo *rdoZDC= (BrZdcRdo*) inNode->GetObject ("RdoZDC");
        if (!rdoZDC) return;

        if (!strcmp (GetName (), "ZDCLeft"))
        {
            bbTime= ((BrBbVertex*) bb)->GetLeftTime ();

            adc [0]= rdoZDC->GetLeftAdc1Lo ();
            adc [1]= rdoZDC->GetLeftAdc2Lo ();
            adc [2]= rdoZDC->GetLeftAdc3Lo ();
            adc [3]= rdoZDC->GetLeftAdcSumLo ();

            tdc [0]= rdoZDC->GetLeftTdc1 ();
            tdc [1]= rdoZDC->GetLeftTdc2 ();
            tdc [2]= rdoZDC->GetLeftTdc3 ();
            tdc [3]= rdoZDC->GetLeftTdcSum ();

            tdcGain [0]= fCalibration->GetTdcGain (1); // 0.100335;
            tdcGain [1]= fCalibration->GetTdcGain (2); // 0.100422;
            tdcGain [2]= fCalibration->GetTdcGain (3); // 0.100524;
            tdcGain [3]= fCalibration->GetTdcGain (4); // 0.099347;

            //fTime [0]->Fill (rdoZDC->GetLeftTime1 ());
            //fTime [1]->Fill (rdoZDC->GetLeftTime2 ());
            //fTime [2]->Fill (rdoZDC->GetLeftTime3 ());
        }

        if (!strcmp (GetName (), "ZDCRight"))
        {
            bbTime= ((BrBbVertex*) bb)->GetRightTime ();

            adc [0]= rdoZDC->GetRightAdc1Lo ();
            adc [1]= rdoZDC->GetRightAdc2Lo ();
            adc [2]= rdoZDC->GetRightAdc3Lo ();
            adc [3]= rdoZDC->GetRightAdcSumLo ();

            tdc [0]= rdoZDC->GetRightTdc1 ();
            tdc [1]= rdoZDC->GetRightTdc2 ();
            tdc [2]= rdoZDC->GetRightTdc3 ();
            tdc [3]= rdoZDC->GetRightTdcSum ();

            tdcGain [0]= fCalibration->GetTdcGain (1); // 0.099233;
            tdcGain [1]= fCalibration->GetTdcGain (2); // 0.099908;
            tdcGain [2]= fCalibration->GetTdcGain (3); // 0.099853;
            tdcGain [3]= fCalibration->GetTdcGain (4); // 0.099986;

            //fTime [0]->Fill (rdoZDC->GetRightTime1 ());
            //fTime [1]->Fill (rdoZDC->GetRightTime2 ());
            //fTime [2]->Fill (rdoZDC->GetRightTime3 ());
        }

        //----- fill slewing profiles
        //----- should I subract pedestals from ADC values in sqrt?

        if (bbTime > 0.0)
        {
            if (adc [0] > 0 && adc [0] < 2000 && tdc [0] > 0 && tdc [0] < 2000) fSlewProf [0]->Fill (1.0/ sqrt (adc [0]), tdc [0]* tdcGain [0]- bbTime);
            if (adc [1] > 0 && adc [1] < 2000 && tdc [1] > 0 && tdc [1] < 2000) fSlewProf [1]->Fill (1.0/ sqrt (adc [1]), tdc [1]* tdcGain [1]- bbTime);
            if (adc [2] > 0 && adc [2] < 2000 && tdc [2] > 0 && tdc [2] < 2000) fSlewProf [2]->Fill (1.0/ sqrt (adc [2]), tdc [2]* tdcGain [2]- bbTime);
            if (adc [3] > 0 && adc [3] < 2000 && tdc [3] > 0 && tdc [3] < 2000) fSlewProf [3]->Fill (1.0/ sqrt (adc [3]), tdc [3]* tdcGain [3]- bbTime);
        }


	if(fUseNtuple){
	  Float_t ntdata[17];

	  ntdata [0]= ((BrBbVertex*) bb)->GetLeftTime ();
	  ntdata [1]= ((BrBbVertex*) bb)->GetRightTime ();
	  ntdata [2]= ((BrBbVertex*) bb)->GetZ0 ();
	  ntdata [3]= ((BrBbVertex*) bb)->GetTime0 ();
	  ntdata [4]= ((BrBbVertex*) bb)->GetMethod ();
	  
	  ntdata [5]= rdoZDC->GetLeftTime1 ();
	  ntdata [6]= rdoZDC->GetLeftTime2 ();
	  ntdata [7]= rdoZDC->GetLeftTime3 ();
	  ntdata [8]= rdoZDC->GetRightTime1 ();
	  ntdata [9]= rdoZDC->GetRightTime2 ();
	  ntdata [10]= rdoZDC->GetRightTime3 ();
	  
	  ntdata [11]= rdoZDC->GetLeftAdc1Lo ();
	  ntdata [12]= rdoZDC->GetLeftAdc2Lo ();
	  ntdata [13]= rdoZDC->GetLeftAdc3Lo ();
	  ntdata [14]= rdoZDC->GetRightAdc1Lo ();
	  ntdata [15]= rdoZDC->GetRightAdc2Lo ();
	  ntdata [16]= rdoZDC->GetRightAdc3Lo ();
	  
	  fNtuple->Fill (ntdata);
	}
    }
}



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



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

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

    Int_t i, j;

    Float_t min1, min2, min3, minSum;		// fitting range:
    Float_t max1, max2, max3, maxSum;	// to be determined phenomenologically
    Float_t slewPar [5][4];			// 5 slew parameters for 4 modules (1, 2, 3, sum)

    TF1 *slewFunc [4];

    //----- fitting ranges here were specifically adjusted up to set of runs
    //----- run005677seq000.dat - run005685seq007.dat from /brahms/data05/data/raw

    //----- fit with pol3 function

    if (!strcmp (GetName (), "ZDCLeft"))
    {
        min1= 0.042; max1= 0.15;
        min2= 0.045; max2= 0.135;
        min3= 0.045; max3= 0.127;
        minSum= 0.035; maxSum= 0.12;
    }

    if (!strcmp (GetName (), "ZDCRight"))
    {
        min1= 0.04; max1= 0.14;
        min2= 0.045; max2= 0.15;
        min3= 0.03; max3= 0.13;
        minSum= 0.035; maxSum= 0.125;
    }

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

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

    for (j= 0; j < fNoModules; j++)
    {
        for (i= 0; i < 5; i++) slewPar [i][j]= slewFunc [j]->GetParameter (i);
    }

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

    for (j= 0; j < fNoModules; j++)
    {
        fCalibration->SetSlewpar1 (j+1, slewPar [0][j]);
        fCalibration->SetSlewpar2 (j+1, slewPar [1][j]);
        fCalibration->SetSlewpar3 (j+1, slewPar [2][j]);
        fCalibration->SetSlewpar4 (j+1, slewPar [3][j]);
        fCalibration->SetSlewpar5 (j+1, slewPar [4][j]);
    }

    fCalibration->SetComment ("Slewpar1", "Generated by  BrZdcSlewCalModule: fit with a pol3 function");
    fCalibration->SetComment ("Slewpar2", "Generated by  BrZdcSlewCalModule: fit with a pol3 function");
    fCalibration->SetComment ("Slewpar3", "Generated by  BrZdcSlewCalModule: fit with a pol3 function");
    fCalibration->SetComment ("Slewpar4", "Generated by  BrZdcSlewCalModule: fit with a pol3 function");
    fCalibration->SetComment ("Slewpar5", "Generated by  BrZdcSlewCalModule: fit with a pol3 function");

    if (fSaveAscii) SaveAscii ();
}



 void BrZdcSlewCalModule::SaveAscii ()
{
    Int_t i;
    Float_t par1, par2, par3, par4, par5;

    //----- save parameters to an ASCII file

    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); //----- fixed mode seems to provide better precision than scientific
    
    file << "********************************************************" << endl;
    file << "*" << endl;
    file << "*" << "	" << GetName () << " slew correction: events from run #" << run->GetRunNo () << endl;
    file << "*" << endl;
    file << "********************************************************" << endl;
    file << "*" << endl;
    file << "*" << "	";
    
    file.width (8);
    file << "Mod.#";
    
    file.width (15);
    file << "Slewpar1";
    
    file.width (15);
    file << "Slewpar2";
    
    file.width (15);
    file << "Slewpar3";
    
    file.width (15);
    file << "Slewpar4";
    
    file.width (15);
    file << "Slewpar5" << endl;
    
    file << "*" << endl;
 
    for (i= 0; i < fNoModules; i++)
    {
        par1= fCalibration->GetSlewpar1 (i+1);
        par2= fCalibration->GetSlewpar2 (i+1);
        par3= fCalibration->GetSlewpar3 (i+1);
        par4= fCalibration->GetSlewpar4 (i+1);
        par5= fCalibration->GetSlewpar5 (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;
        
        file.width (15);
        file << par5 << endl;
    }
}



 void BrZdcSlewCalModule::ReadAscii ()
{
    //----- save parameters 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, par5; 
    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 >> par5;

        fCalibration->SetSlewpar1 (i+1, par1);
        fCalibration->SetSlewpar2 (i+1, par2);
        fCalibration->SetSlewpar3 (i+1, par3);
        fCalibration->SetSlewpar4 (i+1, par4);
        fCalibration->SetSlewpar5 (i+1, par5);
    }
}

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