]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDCalibChamberStatus.cxx
TRD pretrigger simulation based on digits of T0, V0 and TOF, simulation control is...
[u/mrichter/AliRoot.git] / TRD / AliTRDCalibChamberStatus.cxx
CommitLineData
6c1053a8 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16////////////////////////////////////////////////////////////////////////////
17// //
18// AliTRDCalibChamberStatus: to determine which half chambers are off //
19// Produce a AliTRDCalChamberStatus calibration object //
20// Check with the AliTRDCalDCSFEE info //
21// //
22// //
23// Authors: //
24// J. Book (jbook@ikf.uni-frankfurt.de) //
25// R. Bailhache (rbailhache@ikf.uni-frankfurt.de) //
26// //
27////////////////////////////////////////////////////////////////////////////
28
29
30//Root includes
31#include <THnSparse.h>
32
33#include <TDirectory.h>
34#include <TFile.h>
35
36//AliRoot includes
37#include "AliRawReader.h"
38
39//header file
40#include "AliLog.h"
41#include "AliTRDCalibChamberStatus.h"
42#include "AliTRDgeometry.h"
43#include "AliTRDdigitsManager.h"
44#include "AliTRDSignalIndex.h"
45#include "AliTRDrawFastStream.h"
46#include "./Cal/AliTRDCalChamberStatus.h"
47#include "./Cal/AliTRDCalDCS.h"
48#include "./Cal/AliTRDCalDCSFEE.h"
49
50#ifdef ALI_DATE
51#include "event.h"
52#endif
53
54ClassImp(AliTRDCalibChamberStatus) /*FOLD00*/
55
56//_____________________________________________________________________
57AliTRDCalibChamberStatus::AliTRDCalibChamberStatus() : /*FOLD00*/
58 TObject(),
59 fDetector(-1),
60 fNumberOfTimeBins(0),
61 fCounterEventNotEmpty(0),
62 fCalChamberStatus(0x0),
63 fHnSparseI(0x0),
64 fHnSparseHCM(0x0),
65 fHnSparseEvtDet(0x0),
66 fHnSparseDebug(0x0),
67 fHnSparseMCM(0x0),
68 fDebugLevel(0)
69{
70 //
71 // default constructor
72 //
73
74}
75//_____________________________________________________________________
76AliTRDCalibChamberStatus::AliTRDCalibChamberStatus(const AliTRDCalibChamberStatus &ped) : /*FOLD00*/
77 TObject(ped),
78 fDetector(ped.fDetector),
79 fNumberOfTimeBins(ped.fNumberOfTimeBins),
80 fCounterEventNotEmpty(ped.fCounterEventNotEmpty),
81 fCalChamberStatus(ped.fCalChamberStatus),
82 fHnSparseI(ped.fHnSparseI),
83 fHnSparseHCM(ped.fHnSparseHCM),
84 fHnSparseEvtDet(ped.fHnSparseEvtDet),
85 fHnSparseDebug(ped.fHnSparseDebug),
86 fHnSparseMCM(ped.fHnSparseMCM),
87 fDebugLevel(ped.fDebugLevel)
88{
89 //
90 // copy constructor
91 //
92
93}
94//_____________________________________________________________________
95AliTRDCalibChamberStatus& AliTRDCalibChamberStatus::operator = (const AliTRDCalibChamberStatus &source)
96{
97 //
98 // assignment operator
99 //
100 if (&source == this) return *this;
101 new (this) AliTRDCalibChamberStatus(source);
102
103 return *this;
104}
105//_____________________________________________________________________
106AliTRDCalibChamberStatus::~AliTRDCalibChamberStatus() /*FOLD00*/
107{
108 //
109 // destructor
110 //
111 if(fCalChamberStatus){
112 delete fCalChamberStatus;
113 }
114 if(fHnSparseI) {
115 delete fHnSparseI;
116 }
117 if(fHnSparseHCM) {
118 delete fHnSparseHCM;
119 }
120 if(fHnSparseEvtDet) {
121 delete fHnSparseEvtDet;
122 }
123 if(fHnSparseDebug) {
124 delete fHnSparseDebug;
125 }
126 if(fHnSparseMCM) {
127 delete fHnSparseMCM;
128 }
129}
130
131//_____________________________________________________________________
132void AliTRDCalibChamberStatus::Init()
133{
134 //
135 // Init the different THnSparse
136 //
137 //
138
139 //
140 // Init the fHnSparseI
141 //
142
143 //create the map
144 Int_t thnDimEvt[4]; // sm, layer, stack, halfchamber
145 thnDimEvt[0] = 18;
146 thnDimEvt[1] = 6;
147 thnDimEvt[2] = 5;
148 thnDimEvt[3] = 2;
149 //arrays for lower bounds :
150 Double_t* binEdgesEvt[4];
151 for(Int_t ivar = 0; ivar < 4; ivar++)
152 binEdgesEvt[ivar] = new Double_t[thnDimEvt[ivar] + 1];
153 //values for bin lower bounds
154 for(Int_t i=0; i<=thnDimEvt[0]; i++) binEdgesEvt[0][i]= 0.0 + (18.0)/thnDimEvt[0]*(Double_t)i;
155 for(Int_t i=0; i<=thnDimEvt[1]; i++) binEdgesEvt[1][i]= 0.0 + (6.0)/thnDimEvt[1]*(Double_t)i;
156 for(Int_t i=0; i<=thnDimEvt[2]; i++) binEdgesEvt[2][i]= 0.0 + (5.0)/thnDimEvt[2]*(Double_t)i;
157 for(Int_t i=0; i<=thnDimEvt[3]; i++) binEdgesEvt[3][i]= 0.0 + (2.0)/thnDimEvt[3]*(Double_t)i;
158
159 //create the THnSparse
160 fHnSparseI = new THnSparseI("NumberOfEntries","NumberOfEntries",4,thnDimEvt);
161 for (int k=0; k<4; k++) {
162 fHnSparseI->SetBinEdges(k,binEdgesEvt[k]);
163 }
164 fHnSparseI->Sumw2();
165
166 //
167 // Init the fHnSparseHCM (THnSparseI)
168 //
169
170 //create the THnSparse
171 fHnSparseHCM = new THnSparseI("HCMerrors","HCMerrors",4,thnDimEvt);
172 for (int k=0; k<4; k++) {
173 fHnSparseHCM->SetBinEdges(k,binEdgesEvt[k]);
174 }
175 fHnSparseHCM->Sumw2();
176
177
178 //---------//
179 // Debug //
180 if(fDebugLevel > 0) {
181
182 //
183 // Init the fHnSparseEvtDet (THnSparseI)
184 //
185
186 //create the map
187 Int_t thnDimEvts[3]; // event, detector, halfchamber
188 thnDimEvts[0] = 10000;
189 thnDimEvts[1] = 540;
190 thnDimEvts[2] = 2;
191 //arrays for lower bounds :
192 Double_t* binEdgesEvts[3];
193 for(Int_t ivar = 0; ivar < 3; ivar++)
194 binEdgesEvts[ivar] = new Double_t[thnDimEvts[ivar] + 1];
195 //values for bin lower bounds
196 for(Int_t i=0; i<=thnDimEvts[0]; i++) binEdgesEvts[0][i]= 0.0 + (10000.0)/thnDimEvts[0]*(Double_t)i;
197 for(Int_t i=0; i<=thnDimEvts[1]; i++) binEdgesEvts[1][i]= 0.0 + (540.0)/thnDimEvts[1]*(Double_t)i;
198 for(Int_t i=0; i<=thnDimEvts[2]; i++) binEdgesEvts[2][i]= 0.0 + (2.0)/thnDimEvts[2]*(Double_t)i;
199
200 //create the THnSparse
201 fHnSparseEvtDet = new THnSparseI("NumberOfEntriesPerEvent","NumberOfEntriesPerEvent",3,thnDimEvts);
202 for (int k=0; k<3; k++) {
203 fHnSparseEvtDet->SetBinEdges(k,binEdgesEvts[k]);
204 }
205 fHnSparseEvtDet->Sumw2();
206
207 //
208 // Init the fHnSparseDebug (THnSparseI)
209 //
210
211 //create the THnSparse
212 fHnSparseDebug = new THnSparseI("NumberOfDifferentDecisions","NumberOfDifferentDecisions",4,thnDimEvt);
213 for (int k=0; k<4; k++) {
214 fHnSparseDebug->SetBinEdges(k,binEdgesEvt[k]);
215 }
216 fHnSparseDebug->Sumw2();
217
218 //
219 // Init the fHnSparseMCM (THnSparseI)
220 //
221
222 //create the map
223 Int_t thnDimEvtt[6]; // sm, layer, stack, ROB, MCM
224 thnDimEvtt[0] = 18;
225 thnDimEvtt[1] = 6;
226 thnDimEvtt[2] = 5;
227 thnDimEvtt[3] = 8;
228 thnDimEvtt[4] = 16;
229 thnDimEvtt[5] = 16;
230 //arrays for lower bounds :
231 Double_t* binEdgesEvtt[6];
232 for(Int_t ivar = 0; ivar < 6; ivar++)
233 binEdgesEvtt[ivar] = new Double_t[thnDimEvtt[ivar] + 1];
234 //values for bin lower bounds
235 for(Int_t i=0; i<=thnDimEvtt[0]; i++) binEdgesEvtt[0][i]= 0.0 + (18.0)/thnDimEvtt[0]*(Double_t)i;
236 for(Int_t i=0; i<=thnDimEvtt[1]; i++) binEdgesEvtt[1][i]= 0.0 + (6.0)/thnDimEvtt[1]*(Double_t)i;
237 for(Int_t i=0; i<=thnDimEvtt[2]; i++) binEdgesEvtt[2][i]= 0.0 + (5.0)/thnDimEvtt[2]*(Double_t)i;
238 for(Int_t i=0; i<=thnDimEvtt[3]; i++) binEdgesEvtt[3][i]= 0.0 + (8.0)/thnDimEvtt[3]*(Double_t)i;
239 for(Int_t i=0; i<=thnDimEvtt[4]; i++) binEdgesEvtt[4][i]= 0.0 + (16.0)/thnDimEvtt[4]*(Double_t)i;
240 for(Int_t i=0; i<=thnDimEvtt[5]; i++) binEdgesEvtt[5][i]= 0.0 + (16.0)/thnDimEvtt[5]*(Double_t)i;
241
242 //create the THnSparse
243 fHnSparseMCM = new THnSparseI("MCMerrorDCS","MCMerrorDCS",6,thnDimEvtt);
244 for (int k=0; k<6; k++) {
245 fHnSparseMCM->SetBinEdges(k,binEdgesEvtt[k]);
246 }
247 fHnSparseMCM->Sumw2();
248
249 }
250 // Debug //
251 //---------//
252
253}
254//_____________________________________________________________________
255void AliTRDCalibChamberStatus::ProcessEvent(AliRawReader * rawReader, Int_t nevents_physics)
256{
257 //
258 // Event Processing loop
259 //
260 //
261
262 Bool_t notEmpty = kFALSE;
263
264 AliTRDrawFastStream *rawStream = new AliTRDrawFastStream(rawReader);
1a278748 265 rawStream->SetNoErrorWarning();
6c1053a8 266 rawStream->SetSharedPadReadout(kFALSE);
267
268 AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager(kTRUE);
269 digitsManager->CreateArrays();
270
271 Int_t det = 0;
272 while ((det = rawStream->NextChamber(digitsManager, NULL, NULL)) >= 0) {
273
274 //nextchamber loop
275
276 // do the QA analysis
277 if (digitsManager->GetIndexes(det)->HasEntry()) {//QA
278 // printf("there is ADC data on this chamber!\n");
279
280 AliTRDSignalIndex *indexes = digitsManager->GetIndexes(det);
281 if (indexes->IsAllocated() == kFALSE) {
282 // AliError("Indexes do not exist!");
283 break;
284 }
285
286 Int_t iRow = 0;
287 Int_t iCol = 0;
288 indexes->ResetCounters();
289
290 while (indexes->NextRCIndex(iRow, iCol)){
291 Int_t iMcm = (Int_t)(iCol/18); // current group of 18 col pads
292
293 Int_t layer = AliTRDgeometry::GetLayer(det);
294 Int_t sm = AliTRDgeometry::GetSector(det);
295 Int_t stac = AliTRDgeometry::GetStack(det);
296 Double_t rphi = 0.5;
297 if(iMcm > 3) rphi = 1.5;
298
299 Double_t val[4] = {sm,layer,stac,rphi};
300 fHnSparseI->Fill(&val[0]);
301 notEmpty = kTRUE;
302
303 //---------//
304 // Debug //
305 if(fDebugLevel > 0) {
306 Int_t detector = AliTRDgeometry::GetDetector(layer,stac,sm);
307 Double_t valu[3] = {nevents_physics,detector,rphi};
308 fHnSparseEvtDet->Fill(&valu[0]);
309 }
310 // Debug //
311 //---------//
312 }
313
314 }
315 digitsManager->ClearArrays(det);
316 }
317
318 if(notEmpty) fCounterEventNotEmpty++;
319
320 if(digitsManager) delete digitsManager;
321 if(rawStream) delete rawStream;
322
323}
324//_____________________________________________________________________
325Bool_t AliTRDCalibChamberStatus::TestEventHisto(Int_t nevent) /*FOLD00*/
326{
327 //
328 // Test event loop
329 // fill the fHnSparseI with entries
330 //
331
332 AliTRDgeometry geo;
333
334
335 for(Int_t ievent=0; ievent<nevent; ievent++){
336 for (Int_t ism=0; ism<18; ism++){
337 for (Int_t istack=0; istack<5; istack++){
338 for (Int_t ipl=0; ipl<6; ipl++){
339 for (Int_t icol=0; icol<geo.GetColMax(ipl); icol++){
340 Int_t side = 0;
341 if(icol > 72) side = 1;
342 Double_t val[4] = {ism,ipl,istack,side};
343 fHnSparseI->Fill(&val[0]);
344 }
345 }
346 }
347 }
348 }
349
350 return kTRUE;
351
352}
353//_____________________________________________________________________
354void AliTRDCalibChamberStatus::AnalyseHisto() /*FOLD00*/
355{
356 //
357 // Create the AliTRDCalChamberStatus according to the fHnSparseI
358 //
359
360 if(fCalChamberStatus) delete fCalChamberStatus;
361 fCalChamberStatus = new AliTRDCalChamberStatus();
362
363 // Check if enough events to say something
364 if(fCounterEventNotEmpty < 30) {
365 // Say all installed
366 for (Int_t ism=0; ism<18; ism++) {
367 for (Int_t ipl=0; ipl<6; ipl++) {
368 for (Int_t istack=0; istack<5; istack++) {
369 // layer, stack, sector
370 Int_t det = AliTRDgeometry::GetDetector(ipl,istack,ism);
371 fCalChamberStatus->SetStatus(det,1);
372 }
373 }
374 }
375 return;
376 }
377
378 // Mask out all chambers
379 for (Int_t ism=0; ism<18; ism++) {
380 for (Int_t ipl=0; ipl<6; ipl++) {
381 for (Int_t istack=0; istack<5; istack++) {
382 // layer, stack, sector
383 Int_t det = AliTRDgeometry::GetDetector(ipl,istack,ism);
384 fCalChamberStatus->SetStatus(det,2);
385 }
386 }
387 }
388
389 // Unmask good chambers
390 Int_t coord[4];
391 for(Int_t bin = 0; bin < fHnSparseI->GetNbins(); bin++) {
392
393 fHnSparseI->GetBinContent(bin,coord);
394 // layer, stack, sector
395 Int_t detector = AliTRDgeometry::GetDetector(coord[1]-1,coord[2]-1,coord[0]-1);
396
397 //
398 // Check which halfchamber side corresponds to the bin number (0=A, 1=B)
399 // Change the status accordingly
400 //
401
402 switch(fCalChamberStatus->GetStatus(detector))
403 {
404 case 1: break; // no changes
405 case 2:
406 if(coord[3]-1==0) {
407 fCalChamberStatus->SetStatus(detector,4); break; // only SideB is masked
408 }
409 else {
410 fCalChamberStatus->SetStatus(detector,3); break; // only SideA is masked
411 }
412 case 3: fCalChamberStatus->SetStatus(detector,1); break; // unmask SideA
413 case 4: fCalChamberStatus->SetStatus(detector,1); break; // unmask SideB
414 }
415 }
416
417
418}
419//_____________________________________________________________________
420void AliTRDCalibChamberStatus::CheckEORStatus(AliTRDCalDCS *calDCS) /*FOLD00*/
421{
422 //
423 // Correct the AliTRDCalChamberStatus according to the AliTRDCalDCS
424 // Using globale state of the HalfChamberMerger (HCM)
425 //
426
427 for(Int_t det = 0; det < 540; det++) {
428 AliTRDCalDCSFEE* calDCSFEEEOR = calDCS->GetCalDCSFEEObj(det);
429 if(!calDCSFEEEOR) { continue;}
430
431 // MCM Global State Machine State Definitions
432 // low_power = 0,
433 // test = 1,
434 // wait_pre = 3,
435 // preproc = 7,
436 // zero_sp = 8,
437 // full_rd = 9,
438 // clear_st = 11,
439 // wait_L1 = 12,
440 // tr_send = 14,
441 // tr_proc = 15
442
443 Int_t sm = AliTRDgeometry::GetSector(det);
444 Int_t lay = AliTRDgeometry::GetLayer(det);
445 Int_t stac = AliTRDgeometry::GetStack(det);
446
447 Int_t stateA = calDCSFEEEOR->GetMCMGlobalState(4,17); // HCM Side A
448 Int_t stateB = calDCSFEEEOR->GetMCMGlobalState(5,17); // HCM Side B
449 Int_t rphi = -1;
450
451 //printf("DCS: stateA %d \t stateB %d \n",stateA,stateB);
452 if(stateA!=3 && stateA!=9) rphi = 1;
453 Double_t vals[4] = {sm,lay,stac,rphi};
454 if(rphi!=-1) fHnSparseHCM->Fill(&vals[0]);
455
456 if(stateB!=3 && stateB!=9) rphi = 2;
457 vals[3] = rphi;
458 if(rphi!=-1) fHnSparseHCM->Fill(&vals[0]);
1a278748 459
6c1053a8 460 //---------//
461 // Debug //
462 if(fDebugLevel > 0) {
a90e31fc 463 if( ((fCalChamberStatus->GetStatus(det) <= 1) && (stateA!=3 && stateA!=9)) ||
464 ((fCalChamberStatus->GetStatus(det) <= 1) && (stateB!=3 && stateB!=9)) ||
1a278748 465 ((fCalChamberStatus->GetStatus(det) == 2) && (stateA==3 || stateA==9)) ||
466 ((fCalChamberStatus->GetStatus(det) == 2) && (stateB==3 || stateB==9)) ||
467 ((fCalChamberStatus->GetStatus(det) == 3) && (stateA==3 || stateA==9)) ||
468 ((fCalChamberStatus->GetStatus(det) == 3) && (stateB!=3 && stateB!=9)) ||
469 ((fCalChamberStatus->GetStatus(det) == 4) && (stateB==3 || stateB==9)) ||
470 ((fCalChamberStatus->GetStatus(det) == 4) && (stateA!=3 && stateA!=9)) )
6c1053a8 471 {
472 //printf(" Different half chamber status in DCS and DATA!!\n");
473 Double_t val[4] = {sm,lay,stac,1};
474 fHnSparseDebug->Fill(&val[0]);
475
476 if(rphi!=-1) { // error in DCS information
477 // Fill MCM status map
478 for(Int_t ii = 0; ii < 8; ii++) { //ROB loop
479 for(Int_t i = 0; i < 16; i++) { //MCM loop
480 Double_t valss[6] = {sm,lay,stac,ii,i
481 ,calDCSFEEEOR->GetMCMGlobalState(ii,i)};
482 fHnSparseMCM->Fill(&valss[0]);
483 }
484 }
485 }
486 }
487 }
488 //---------//
489 // Debug //
1a278748 490
491 //////////////////////////////////////////
492 // Change the status according to DCS
493 ///////////////////////////////////////////
494
495 // First put bad status if seen by DCS
496 /////////////////////////////////////////
497 if(stateA!=3 && stateA!=9 && stateB!=3 && stateB!=9) {
498 // completely masked from DCS
499 fCalChamberStatus->SetStatus(det,2);
500 }
501 if((stateA!=3 && stateA!=9) && (stateB==3 || stateB==9)) {
502 // Only A side masked from DCS
503 if(fCalChamberStatus->GetStatus(det) != 3) {
504 // But not from Data
505 if(fCalChamberStatus->GetStatus(det) == 4) fCalChamberStatus->SetStatus(det,2);
506 if(fCalChamberStatus->GetStatus(det) <= 1) fCalChamberStatus->SetStatus(det,3);
507 }
508 }
509 if((stateA==3 || stateA==9) && (stateB!=3 && stateB!=9)) {
510 // Only B side masked from DCS
511 if(fCalChamberStatus->GetStatus(det) != 4) {
512 // But not from Data
513 if(fCalChamberStatus->GetStatus(det) == 3) fCalChamberStatus->SetStatus(det,2);
514 if(fCalChamberStatus->GetStatus(det) <= 1) fCalChamberStatus->SetStatus(det,4);
515 }
516 }
517
518 // Then release in case of error 3 from DCS
519 /////////////////////////////////////////////
520 if(fCalChamberStatus->GetStatus(det)==2) {
521 if((stateA==3) && (stateB==3)) fCalChamberStatus->SetStatus(det,1);
522 if((stateA==3) && (stateB!=3)) fCalChamberStatus->SetStatus(det,4);
523 if((stateA!=3) && (stateB==3)) fCalChamberStatus->SetStatus(det,3);
524 }
525 if((fCalChamberStatus->GetStatus(det)==3) && (stateA==3)) fCalChamberStatus->SetStatus(det,1);
526 if((fCalChamberStatus->GetStatus(det)==4) && (stateB==3)) fCalChamberStatus->SetStatus(det,1);
527
6c1053a8 528 }
529
530}
531
532//_____________________________________________________________________________________
533void AliTRDCalibChamberStatus::Add(AliTRDCalibChamberStatus *calibChamberStatus) /*FOLD00*/
534{
535 //
536 // Add the THnSparseI of this calibChamberStatus
537 //
538
539 fCounterEventNotEmpty += calibChamberStatus->GetNumberEventNotEmpty();
540
541 THnSparseI *hnSparseI = calibChamberStatus->GetSparseI();
542 if(!hnSparseI) return;
543
544 if(!fHnSparseI) {
545 fHnSparseI = (THnSparseI *) hnSparseI->Clone();
546 }
547 else {
548 fHnSparseI->Add(hnSparseI);
549 }
550
551
552}
553//_____________________________________________________________________
554void AliTRDCalibChamberStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
555{
556 //
557 // Write class to file
558 //
559
560 TString sDir(dir);
561 TString option;
562
563 if ( append )
564 option = "update";
565 else
566 option = "recreate";
567
568 TDirectory *backup = gDirectory;
569 TFile f(filename,option.Data());
570 f.cd();
571 if ( !sDir.IsNull() ){
572 f.mkdir(sDir.Data());
573 f.cd(sDir);
574 }
575 this->Write();
576 f.Close();
577
578 if ( backup ) backup->cd();
579}
580