ROOT logo
// $Id$
// Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007

/**************************************************************************
 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
 * full copyright notice.                                                 *
 **************************************************************************/

#if !defined(__CINT__) || defined(__MAKECINT__)
#include <iostream>

#include <TGeoGlobalMagField.h>
#include <TGeoManager.h>
#include <TMap.h>
#include <TString.h>
#include <TSystem.h>
#include <TStyle.h>
#include <TTree.h>
#include <TRegexp.h>
#include <TEveManager.h>
#include <TEveElement.h>
#include <TEvePointSet.h>
#include <TEveTreeTools.h>

#include <AliAltroRawStreamV3.h>
#include <AliEveEventManager.h>
#include <AliEveTPCData.h>
#include <AliEveTPCSector2D.h>
#include <AliEveTPCSector3D.h>
#include <AliGeomManager.h>
#include <AliMagF.h>
#include <AliRawReader.h>

#include <AliGRPObject.h>
#include <AliCDBManager.h>
#include <AliCDBPath.h>
#include <AliCDBEntry.h>

#include <STEER/STEER/AliGRPManager.h>

#include <TPC/AliTPCcalibDB.h>
#include <TPC/AliTPCTransform.h>
#include <TPC/AliTPCclusterMI.h>
#include <TPC/AliTPCClustersRow.h>

#include <TPC/AliTPCParam.h>
#include <TPC/AliTPCRawStreamV3.h>
#include <TPC/AliTPCRecoParam.h>
#include <TPC/AliTPCReconstructor.h>

#include <HLT/BASE/AliHLTOUT.h>
#include <HLT/BASE/AliHLTSystem.h>
#include <HLT/BASE/AliHLTPluginBase.h>

#include <RAW/AliRawHLTManager.h>
#endif

// read TPC clusters compressed by HLT from RawReader and return output in OutputTree
TTree* readHLTClusters(AliRawReader *RawReader);

// add the HLT clusters from clustersTree to visualization
void renderHLTClusters(TTree* clustersTree);


// Macro to visualise rootified raw-data from TPC.
//
// Use tpc_raw(Int_t mode) in order to run it
// Needs that alieve_init() is already called
// mode = 1 - show only 2D sectors
// mode = 2 - show only 3D sectors
// mode = 3 - show both 2D and 3D sectors
void tpc_raw(Int_t mode = 3)
{
  gStyle->SetPalette(1, 0);

  AliEveEventManager::AssertGeometry();
  AliEveEventManager::AssertMagField();
  
  AliRawReader *reader = AliEveEventManager::AssertRawReader();
  reader->Reset();
  
  AliTPCRawStreamV3 input(reader);
  reader->Select("TPC"); // ("TPC", firstRCU, lastRCU);

  AliEveTPCData *x = new AliEveTPCData;
  // x->SetLoadPedestal(5);
  x->SetLoadThreshold(5);
  x->SetAutoPedestal(kTRUE);

  x->LoadRaw(input, kTRUE, kTRUE);

  gEve->DisableRedraw();

  TEveElementList* sec2d = new TEveElementList("TPC 2D");
  gEve->AddElement(sec2d);

  TEveElementList* sec3d = new TEveElementList("TPC 3D");
  gEve->AddElement(sec3d);

  AliEveTPCSector2D *s;
  AliEveTPCSector3D *t;
  
  for (Int_t i=0; i<=35; ++i) {
    if (mode & 1) {
      s = new AliEveTPCSector2D(Form("2D sector %d",i));
      s->SetSectorID(i);
      s->SetAutoTrans(kTRUE); // place on proper 3D coordinates
      s->SetDataSource(x);
      s->SetFrameColor(36);
      sec2d->AddElement(s);
      s->IncRTS();
    }
    if (mode & 2) {
      t = new AliEveTPCSector3D(Form("3D sector %d",i));
      t->SetSectorID(i);
      t->SetAutoTrans(kTRUE);
      t->SetDataSource(x);
      sec3d->AddElement(t);
      t->IncRTS();
    }
  }
  
  // Display TPC clusters compressed by HLT
  TTree* hltClustersTree = readHLTClusters(reader); // read HLT compressed clusters from TPC from raw reader and output them in hltClustersTree
  if(hltClustersTree) renderHLTClusters(hltClustersTree);
  
  gEve->EnableRedraw();
  gEve->Redraw3D();
  
}

// read TPC Clusters compressed by HLT from a raw filename
TTree* readHLTClusters(AliRawReader *fRawReader)
{
/*************************
*  Raw IO Stuff
**************************/
  AliRawReader* rawReader = AliRawHLTManager::CreateRawReaderHLT(fRawReader, "TPC");
  rawReader->Select("TPC");

/*************************
*  HLT IO Stuff
**************************/
  AliHLTSystem* hltSystem=AliHLTPluginBase::GetInstance();

  AliHLTOUT* hltOut = AliHLTOUT::New(rawReader); 
  hltOut->Init();
  hltSystem->InitHLTOUT(hltOut);
  
/*************************
*  GRP Stuff
**************************/
  AliGRPManager manGRP;
  AliRecoParam* recoParam = AliEveEventManager::AssertRecoParams();
  
  AliEveEventManager* curEventMan = (AliEveEventManager*)gEve->GetCurrentEvent();
  
  const AliRunInfo* runInfo = manGRP.GetRunInfo();
  const THashTable* cosmicTriggers = manGRP.GetCosmicTriggers();
  const AliEventInfo* eventinfo = curEventMan->GetEventInfo();
  
  recoParam->SetEventSpecie(runInfo,*eventinfo,cosmicTriggers);
  
  AliTPCRecoParam *tpcRecoParam = (AliTPCRecoParam*)recoParam->GetDetRecoParam(1); // TPC has index 1
  tpcRecoParam->SetUseHLTClusters(3); // reconstruct only HLT clusters
	
/**************************
*  Reconstruction of Clusters
**************************/
  TTree* outputClustersTree = new TTree;
  
  AliTPCReconstructor *reconstructor = new AliTPCReconstructor();
  reconstructor->SetOption("useHLT");
  reconstructor->CreateTracker(); // this will set the option to reconstruct for only HLT clusters
  
  reconstructor->SetRecoParam(tpcRecoParam);
  reconstructor->SetRunInfo((AliRunInfo*)runInfo);
  reconstructor->SetEventInfo((AliEventInfo*)eventinfo);
  
  reconstructor->Reconstruct(rawReader, outputClustersTree);
  
  delete reconstructor;
  
  hltSystem->ReleaseHLTOUT(hltOut);
  
  return outputClustersTree;
}
 
void renderHLTClusters(TTree* clustersTree)
{
 
/**************************
*  Visualization of Clusters
**************************/
  const Int_t kMaxCl=100*160;
  Int_t fNColorBins = 5;
  
  TEvePointSet* clusters = new TEvePointSet(kMaxCl);
  clusters->SetOwnIds(kTRUE);

  TEvePointSetArray * cc = new TEvePointSetArray("TPC Clusters Colorized");
  cc->SetMainColor(kRed);
  cc->SetMarkerStyle(4);
  cc->SetMarkerSize(0.4);
  cc->InitBins("Cluster Charge", fNColorBins, 0., fNColorBins*60.);
  
  cc->GetBin(0)->SetMainColor(kGray);
  cc->GetBin(0)->SetMarkerSize(0.4);
  cc->GetBin(1)->SetMainColor(kBlue);
  cc->GetBin(1)->SetMarkerSize(0.42);
  cc->GetBin(2)->SetMainColor(kCyan);
  cc->GetBin(2)->SetMarkerSize(0.44);
  cc->GetBin(3)->SetMainColor(kGreen);
  cc->GetBin(3)->SetMarkerSize(0.46);
  cc->GetBin(4)->SetMainColor(kYellow);
  cc->GetBin(4)->SetMarkerSize(0.48);
  cc->GetBin(5)->SetMainColor(kRed);
  cc->GetBin(5)->SetMarkerSize(0.50);
  cc->GetBin(6)->SetMainColor(kMagenta);
  cc->GetBin(6)->SetMarkerSize(0.52);
  

// Loop over clusters
  Int_t nentries = clustersTree->GetEntriesFast();
	
  AliTPCClustersRow *clrow = new AliTPCClustersRow();
  clrow->SetClass("AliTPCclusterMI");
  //clrow->SetArray(kMaxCl);
  clustersTree->SetBranchAddress("Segment", &clrow);
  
  for (Int_t i=0; i<nentries; i++) {
    if (!clustersTree->GetEvent(i)) continue;
    
    TClonesArray *cl = clrow->GetArray();
    Int_t ncl = cl->GetEntriesFast();

    while (ncl--)
    {
      AliTPCclusterMI* clusterMI = (AliTPCclusterMI*) cl->At(ncl);

      AliCluster *c = (AliCluster*) cl->UncheckedAt(ncl);
      Float_t g[3]; //global coordinates
      c->GetGlobalXYZ(g);
   
      cc->Fill(g[0], g[1], g[2], clusterMI->GetQ());
      clusters->SetNextPoint(g[0], g[1], g[2]);
      AliCluster *atp = new AliCluster(*clusterMI);
      clusters->SetPointId(atp);

    }
    cl->Clear();
  }

  delete clrow;
	  
  clusters->SetName("TPC Clusters");
  clusters->SetTitle(Form("N=%d", clusters->Size()));

  const TString viz_tag("REC Clusters TPC"); // to be changed
  clusters->ApplyVizTag(viz_tag, "Clusters");
  
  cc->SetRnrSelf(kTRUE);

  gEve->AddElement(cc);

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