BRAT 2.4.5
Class index
Full class index
brahmlib
BRAHMS
ROOT page
// $Id: BrTpcCluster.cxx,v 1.7 2002/06/07 16:09:23 videbaek Exp $
// $Author: videbaek $
// $Date: 2002/06/07 16:09:23 $
//

//_____________________________________________________________________
//
// BrTpcCluster is a BRAHMS data class for storing information for 
// cluster island. This is  used in the pattern 
// recognition phase of local tracking (see BrTPCLocalTrackingModule ).
//
// - references a list of BrTpcSequence that makes up an island.
// - is the higher level object (logically) in tracking of the BrTPCHitCluster.
// The Class is based on MJs Cluster Class.
// 
//
#ifndef WIN32
#include <iostream>
#include <iomanip>
#else
#include <iostream.h>
#include <iomanip.h>
#endif

#include "BrTpcCluster.h"

#ifndef BRAT_BrTpcSequence
#include "BrTpcSequence.h"
#endif

#ifndef ROOT_TList
#include "TList.h"
#endif


//_____________________________________________________________________

ClassImp(BrTpcCluster)


//_____________________________________________________________________
   Int_t BrTpcCluster::Compare(const TObject *clus_o) const
{
//This routine overloads TObject::Compare(TObject *object)
//For BRAHMS, this routine is typically called by QSort
//It needs to return 1 if you want clus_o to be earlier
//in the list and -1 if you want clus_o to be later in
//the list.
//Everything is equal!!!
return 0;
}

//_____________________________________________________________________
   Int_t BrTpcCluster::Compare( TObject *clus_o) 
{
//This routine overloads TObject::Compare(TObject *object)
//For BRAHMS, this routine is typically called by QSort
//It needs to return 1 if you want clus_o to be earlier
//in the list and -1 if you want clus_o to be later in
//the list.
//Everything is equal!!!
return 0;
}

//_____________________________________________________________________
 BrTpcCluster::BrTpcCluster(){
  fNpads     = 0;
  fFirstPad  = 0;
  fNbins     = 0;
  fFirstBin  = 0;
  fNseqs     = 0;
  fMaxADC    = 0;
  fBuild     = kFALSE;
  fSeqList   = 0;
}

//_____________________________________________________________________
 BrTpcCluster::BrTpcCluster( Int_t rowNum, BrTpcSequence &firstSeq )
  :fRow(rowNum)
{
  fNpads     = 0;
  fFirstPad  = 0;
  fNbins     = 0;
  fFirstBin  = 0;
  fNseqs     = 0;
  fMaxADC    = 0;
  fBuild     = kFALSE;
  fSeqList   = 0;
  AddSeq(firstSeq);

}
//____________________________________________________________________
 BrTpcCluster::~BrTpcCluster()
{
  if( fSeqList ) {
    fSeqList->Delete();
    delete fSeqList;
  }
}

//____________________________________________________________________
 BrTpcSequence *BrTpcCluster::GetSeq( Int_t seqNum ) 
{
  BrTpcSequence *seq = (BrTpcSequence *)fSeqList->At( seqNum );
  return seq;
}
//____________________________________________________________________
 void BrTpcCluster::AddSeq(BrTpcSequence &seq, Bool_t copy)
{
  // Add a new sequence to the TList of BrTpcSequences 
  // that make up this cluster. At present done by copy.
  // Should be done by reference only. The TList should not own
  // the Sequence data.
  // if copy = kFALSE don't copy just add (only used when deconvoluting)
  // default is to copy
  if(!fSeqList)
    fSeqList = new TList();

  if(copy) {

    BrTpcSequence *copy = new BrTpcSequence(seq);
    fSeqList->Add(copy);
  } else {
    
    fSeqList->Add(&seq);
  }
  
  fNseqs++;
  fBuild = kFALSE;
}

//____________________________________________________________________
 void BrTpcCluster::Print(Option_t* option) const 
{
  cout << "TPC Cluster id = " << fID 
       << (fBuild ? " " : " not ") << "build" << endl
       << " RowNo     = " << setw(4) << fRow
       << " Npads     = " << setw(4) << fNpads
       << " Sequences = " << setw(4) << fNseqs << endl;
  if(fBuild) 
    cout   << " Mean pad  = " << setw(8) << fMeanPad 
	   << " Var pad   = " << setw(8) << fVarPad << endl 
	   << " Mean time = " << setw(8) << fMeanTime
	   << " Var time  = " << setw(8) << fVarTime << endl
	   << " Adc max   = " << setw(8) << fMaxADC 
	   << " Adc sum   = " << setw(8) << fEnergy << endl;
}
//____________________________________________________________________
 void BrTpcCluster::Build( void )
{
  // This routine fills all the derived quantities for a cluster
  Int_t minPad = 2000, maxPad = 0, minTime = 2000, maxTime = 0; 
  Int_t maxADC = -10; 
  // Has to Double_t unless we want to scale cause the thing breaks
  // for a few clusters otherwise => the old DTime = 0.0.
  Double_t  sum          = 0.0;
  Double_t  pav          = 0.0;
  Double_t  tav          = 0.0;
  Double_t  psig         = 0.0;
  Double_t  tsig         = 0.0;
  
  Int_t i;
  for( i = 0; i < fNseqs; i++ ) {
    BrTpcSequence * const seq = GetSeq( i );
    Int_t nSeq = seq->GetNseq();
    Int_t pad = seq->GetPad();
    Int_t timeStart = seq->GetTime();
    Int_t timeStop = timeStart + seq->GetNseq() - 1;
    Short_t *adcSh = seq->GetSequence();
    if( pad < minPad )
      minPad = pad;
    if( pad > maxPad )
      maxPad = pad;
    if( timeStart < minTime )
      minTime = timeStart;
    if( timeStop > maxTime )
      maxTime = timeStop;
    Int_t j;
    for( j = 0; j < nSeq; j++ ) {
      Double_t adc = adcSh[j];
      Double_t time = timeStart+j;
      sum    +=adc;
      pav    +=adc*pad;
      tav    +=adc*time;
      psig   +=adc*pad*pad;
      tsig   +=adc*time*time;
      if( adc > maxADC )
	maxADC = adcSh[j];
    }
  }
  pav  = pav/sum;
  tav  = tav/sum;
  psig = psig/sum - pav*pav;
  if( psig > 0 )
    psig = TMath::Sqrt( psig );
  tsig = tsig/sum - tav*tav;
  if( tsig > 0 )
    tsig = TMath::Sqrt( tsig );

  if (minPad > pav || pav > maxPad)
    Warning("Build", "Pad mean problem - should not happen!!!");
  if (minTime > tav || tav > maxTime)
    Warning("Build", "Time mean problem - should not happen!!!");
  if (pav-minPad < psig && maxPad-pav < psig) 
    Warning("Build", "Pad width problem - should not happen!!!");
  if (tav-minTime < tsig && maxTime-tav < tsig)  
    Warning("Build", "Time width problem - should not happen!!!");
  if(maxADC < 1 || sum < 1) 
    Warning("Build", "Max ADC or sum is 0 or less - should not happen!!!");

  fMeanPad = pav; 
  fVarPad = psig;
  fMeanTime = tav; 
  fVarTime = tsig;
  fEnergy = (Int_t)sum;
  fFirstPad = minPad;
  fFirstBin = minTime;
  fNpads = maxPad-minPad+1; 
  fNbins = maxTime-minTime+1;
  fMaxADC = maxADC;
  // And very importantly
  fBuild = kTRUE;
}

//____________________________________________________________________
//  $Log: BrTpcCluster.cxx,v $
//  Revision 1.7  2002/06/07 16:09:23  videbaek
//  Print only cluster information if it has been build.
//
//  Revision 1.6  2001/08/15 13:38:49  pchristi
//  Fixed memory leak because of TList being created in Default constructer,
//  so that when clusters were read in for files an extra TList was created for
//  each cluster
//
//  Revision 1.5  2001/07/10 11:19:06  cholm
//  Whoops, small mistake. Fixed.
//
//  Revision 1.4  2001/07/10 11:00:55  cholm
//  Changed the BrTpcCluster::Print method so that it's truely const.  The
//  non-const methods GetNpads() and so on, was called.  This is not
//  allowed by ANSI C++.  To make sure that the user is aware of wether the
//  cluster is build or not, additional output was added.  Changed member
//  kBuild to fBuild.  This was a direct violation of the coding standard
//  in BRAT/ROOT, since kBuild is not a constant.
//
//  Revision 1.3  2001/07/06 10:16:15  pchristi
//  Removed asserts and made warnings instead. Changed addseq so it does not always  have to copy the sequence.
//
//  Revision 1.2  2001/06/25 14:30:36  cholm
//  Made Print conform to TObject
//
//  Revision 1.1  2001/06/22 17:36:22  cholm
//  Change names so that every data class has the same format so that
//  we have that down from the very beginning, so that we will not have to
//  worry about that later on.  The affected classes are:
//
//
//          BrRdoBB             ->        BrBbRdo
//          BrRdoZDC            ->        BrZdcRdo
//  	BrTPCCluster	    ->	      BrTpcCluster
//  	BrTPCClusterTable   ->	      BrTpcClusterTable
//
//  Revision 1.1.1.1  2001/06/21 14:55:03  hagel
//  Initial revision of brat2
//
//  Revision 1.9  2001/06/17 17:20:34  pchristi
//  Added print method, comments, and changed fBuild to non-persistent in cluster.
//  GetAdc() in sequence now returns Short_t.
//
//  Revision 1.8  2001/04/14 22:14:15  videbaek
//  added const Compare method
//
//  Revision 1.7  2000/10/02 17:07:58  pchristi
//  Cleaned up headers and source files in the tpc dir
//
//  Revision 1.6  2000/10/02 12:35:17  pchristi
//  Update of the clustering algorithm and related classes.
//
//  Revision 1.5  2000/08/28 18:05:55  videbaek
//  Move CVS log to end of file
//
//  Revision 1.4  2000/03/08 20:18:01  videbaek
//  Change clustering algorithms by adding a new method. BrTpcCluster is redefined.
//  Tracking reorganized with updates by Alv, Peter and Flemming
//
//  Revision 1.3  1999/12/18 21:37:42  videbaek
//  Obsolete for a long timeBrDigTTPC.cxx
//
//  Revision 1.2  1999/02/25 17:26:15  videbaek
//  Add CVS Id and Log flags. Removed commented code
//
//

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