]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibDB.cxx
Adding QA makers for digits (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibDB.cxx
CommitLineData
c5bbaa2c 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// //
19// Class providing the calibration parameters by accessing the CDB //
20// //
21// Request an instance with AliTPCcalibDB::Instance() //
22// If a new event is processed set the event number with SetRun //
1ac191a6 23// Then request the calibration data ////
f5344549 24//
25//
1ac191a6 26// Calibration data:
8cd9634d 27// 0.) Altro mapping
28// Simulation - not yet
29// Reconstruction - AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
30//
1ac191a6 31// 1.) pad by pad calibration - AliTPCCalPad
f5344549 32//
1ac191a6 33// a.) fPadGainFactor
34// Simulation: AliTPCDigitizer::ExecFast - Multiply by gain
35// Reconstruction : AliTPCclustererMI::Digits2Clusters - Divide by gain
f5344549 36//
1ac191a6 37// b.) fPadNoise -
38// Simulation: AliTPCDigitizer::ExecFast
39// Reconstruction: AliTPCclustererMI::FindClusters(AliTPCCalROC * noiseROC)
8cd9634d 40// Noise depending cut on clusters charge (n sigma)
f5344549 41// c.) fPedestal:
42// Simulation: Not used yet - To be impleneted - Rounding to the nearest integer
43// Reconstruction: Used in AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
44// if data taken without zero suppression
45// Currently switch in fRecoParam->GetCalcPedestal();
46//
47// d.) fPadTime0
48// Simulation: applied in the AliTPC::MakeSector - adding offset
49// Reconstruction: AliTPCTransform::Transform() - remove offset
50// AliTPCTransform::Transform() - to be called
51// in AliTPCtracker::Transform()
8cd9634d 52//
53//
54// 2.) Space points transformation:
55//
56// a.) General coordinate tranformation - AliTPCtransform (see $ALICE_ROOT/TPC/AliTPCtransform.cxx)
57// Created on fly - use the other calibration components
58// Unisochronity - (substract time0 - pad by pad)
59// Drift velocity - Currently common drift velocity - functionality of AliTPCParam
60// ExB effect
61// Simulation - Not used directly (the effects are applied one by one (see AliTPC::MakeSector)
62// Reconstruction -
63// AliTPCclustererMI::AddCluster
64// AliTPCtrackerMI::Transform
65// b.) ExB effect calibration -
66// classes (base class AliTPCExB, implementation- AliTPCExBExact.h AliTPCExBFirst.h)
67// a.a) Simulation: applied in the AliTPC::MakeSector -
68// calib->GetExB()->CorrectInverse(dxyz0,dxyz1);
69// a.b) Reconstruction -
70//
71// in AliTPCtransform::Correct() - called calib->GetExB()->Correct(dxyz0,dxyz1)
72//
96305e49 73// 3.) cluster error, shape and Q parameterization
74//
75//
8cd9634d 76//
c5bbaa2c 77///////////////////////////////////////////////////////////////////////////////
78
418bbcaf 79#include <iostream>
80#include <fstream>
81
c5bbaa2c 82
83#include <AliCDBManager.h>
c5bbaa2c 84#include <AliCDBEntry.h>
85#include <AliLog.h>
86
87#include "AliTPCcalibDB.h"
d6834f5f 88#include "AliTPCAltroMapping.h"
418bbcaf 89#include "AliTPCExB.h"
c5bbaa2c 90
91#include "AliTPCCalROC.h"
92#include "AliTPCCalPad.h"
54472e4f 93#include "AliTPCSensorTempArray.h"
418bbcaf 94#include "AliTPCTransform.h"
d6834f5f 95
418bbcaf 96class AliCDBStorage;
97class AliTPCCalDet;
86df2b3a 98//
99//
100
86df2b3a 101#include "TFile.h"
102#include "TKey.h"
103
104#include "TObjArray.h"
105#include "TObjString.h"
106#include "TString.h"
107#include "AliTPCCalPad.h"
0fe7645c 108#include "AliTPCCalibPulser.h"
86df2b3a 109#include "AliTPCCalibPedestal.h"
110#include "AliTPCCalibCE.h"
111
112
113
114
c5bbaa2c 115
116ClassImp(AliTPCcalibDB)
117
118AliTPCcalibDB* AliTPCcalibDB::fgInstance = 0;
119Bool_t AliTPCcalibDB::fgTerminated = kFALSE;
120
121
122//_ singleton implementation __________________________________________________
123AliTPCcalibDB* AliTPCcalibDB::Instance()
124{
125 //
126 // Singleton implementation
127 // Returns an instance of this class, it is created if neccessary
128 //
129
130 if (fgTerminated != kFALSE)
131 return 0;
132
133 if (fgInstance == 0)
134 fgInstance = new AliTPCcalibDB();
135
136 return fgInstance;
137}
138
139void AliTPCcalibDB::Terminate()
140{
141 //
142 // Singleton implementation
143 // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
144 // This function can be called several times.
145 //
146
147 fgTerminated = kTRUE;
148
149 if (fgInstance != 0)
150 {
151 delete fgInstance;
152 fgInstance = 0;
153 }
154}
155
156//_____________________________________________________________________________
e4dce695 157AliTPCcalibDB::AliTPCcalibDB():
158 fRun(-1),
f5344549 159 fTransform(0),
481f877b 160 fExB(0),
e4dce695 161 fPadGainFactor(0),
162 fPadTime0(0),
e4dce695 163 fPadNoise(0),
164 fPedestals(0),
165 fTemperature(0),
d6834f5f 166 fMapping(0),
96305e49 167 fParam(0),
168 fClusterParam(0)
c5bbaa2c 169{
170 //
171 // constructor
172 //
54472e4f 173 //
c5bbaa2c 174 Update(); // temporary
175}
176
177//_____________________________________________________________________________
178AliTPCcalibDB::~AliTPCcalibDB()
179{
180 //
181 // destructor
182 //
68751c2c 183
184 // don't delete anything, CDB cache is active!
185 //if (fPadGainFactor) delete fPadGainFactor;
186 //if (fPadTime0) delete fPadTime0;
68751c2c 187 //if (fPadNoise) delete fPadNoise;
c5bbaa2c 188}
189
190
191//_____________________________________________________________________________
192AliCDBEntry* AliTPCcalibDB::GetCDBEntry(const char* cdbPath)
193{
194 //
195 // Retrieves an entry with path <cdbPath> from the CDB.
196 //
197 char chinfo[1000];
198
68751c2c 199 AliCDBEntry* entry = AliCDBManager::Instance()->Get(cdbPath, fRun);
c5bbaa2c 200 if (!entry)
201 {
202 sprintf(chinfo,"AliTPCcalibDB: Failed to get entry:\t%s ", cdbPath);
203 AliError(chinfo);
204 return 0;
205 }
206 return entry;
207}
208
209
210//_____________________________________________________________________________
211void AliTPCcalibDB::SetRun(Long64_t run)
212{
213 //
214 // Sets current run number. Calibration data is read from the corresponding file.
215 //
216 if (fRun == run)
217 return;
218 fRun = run;
219 Update();
220}
221
222
223
224void AliTPCcalibDB::Update(){
225 //
226 AliCDBEntry * entry=0;
68751c2c 227
228 Bool_t cdbCache = AliCDBManager::Instance()->GetCacheFlag(); // save cache status
229 AliCDBManager::Instance()->SetCacheFlag(kTRUE); // activate CDB cache
230
c5bbaa2c 231 //
232 entry = GetCDBEntry("TPC/Calib/PadGainFactor");
233 if (entry){
68751c2c 234 //if (fPadGainFactor) delete fPadGainFactor;
c5bbaa2c 235 entry->SetOwner(kTRUE);
236 fPadGainFactor = (AliTPCCalPad*)entry->GetObject();
237 }
238 //
239 entry = GetCDBEntry("TPC/Calib/PadTime0");
240 if (entry){
68751c2c 241 //if (fPadTime0) delete fPadTime0;
c5bbaa2c 242 entry->SetOwner(kTRUE);
243 fPadTime0 = (AliTPCCalPad*)entry->GetObject();
244 }
245 //
c5bbaa2c 246 //
247 entry = GetCDBEntry("TPC/Calib/PadNoise");
248 if (entry){
68751c2c 249 //if (fPadNoise) delete fPadNoise;
c5bbaa2c 250 entry->SetOwner(kTRUE);
251 fPadNoise = (AliTPCCalPad*)entry->GetObject();
252 }
8477f500 253
254 entry = GetCDBEntry("TPC/Calib/Pedestals");
255 if (entry){
256 //if (fPedestals) delete fPedestals;
257 entry->SetOwner(kTRUE);
258 fPedestals = (AliTPCCalPad*)entry->GetObject();
259 }
260
54472e4f 261 entry = GetCDBEntry("TPC/Calib/Temperature");
262 if (entry){
263 //if (fTemperature) delete fTemperature;
264 entry->SetOwner(kTRUE);
265 fTemperature = (AliTPCSensorTempArray*)entry->GetObject();
266 }
267
8477f500 268 entry = GetCDBEntry("TPC/Calib/Parameters");
269 if (entry){
54472e4f 270 //if (fPadNoise) delete fPadNoise;
8477f500 271 entry->SetOwner(kTRUE);
a778f7e3 272 fParam = (AliTPCParam*)(entry->GetObject()->Clone());
8477f500 273 }
274
96305e49 275 entry = GetCDBEntry("TPC/Calib/ClusterParam");
276 if (entry){
277 //if (fPadNoise) delete fPadNoise;
278 entry->SetOwner(kTRUE);
279 fClusterParam = (AliTPCClusterParam*)(entry->GetObject()->Clone());
280 }
281
d6834f5f 282 entry = GetCDBEntry("TPC/Calib/Mapping");
283 if (entry){
284 //if (fPadNoise) delete fPadNoise;
285 entry->SetOwner(kTRUE);
286 TObjArray * array = dynamic_cast<TObjArray*>(entry->GetObject());
287 if (array && array->GetEntriesFast()==6){
288 fMapping = new AliTPCAltroMapping*[6];
289 for (Int_t i=0; i<6; i++){
290 fMapping[i] = dynamic_cast<AliTPCAltroMapping*>(array->At(i));
291 }
292 }
293 }
294
295
296
481f877b 297 entry = GetCDBEntry("TPC/Calib/ExB");
298 if (entry) {
299 entry->SetOwner(kTRUE);
300 fExB=dynamic_cast<AliTPCExB*>(entry->GetObject()->Clone());
301 }
302
f5344549 303 if (!fTransform) {
304 fTransform=new AliTPCTransform();
305 }
8477f500 306
c5bbaa2c 307 //
68751c2c 308 AliCDBManager::Instance()->SetCacheFlag(cdbCache); // reset original CDB cache
309
c5bbaa2c 310}
e4dce695 311AliTPCcalibDB::AliTPCcalibDB(const AliTPCcalibDB& org)
312{
313 //
314 // Copy constructor invalid -- singleton implementation
315 //
316 Error("copy constructor","invalid -- singleton implementation");
317}
318
319AliTPCcalibDB& AliTPCcalibDB::operator= (const AliTPCcalibDB& rhs)
320{
321//
322// Singleton implementation - no assignment operator
323//
324 Error("operator =", "assignment operator not implemented");
325 return *this;
326}
327
86df2b3a 328
329
330void AliTPCcalibDB::CreateObjectList(const Char_t *filename, TObjArray *calibObjects)
331{
418bbcaf 332//
333// Create calibration objects and read contents from OCDB
334//
86df2b3a 335 if ( calibObjects == 0x0 ) return;
336 ifstream in;
337 in.open(filename);
338 if ( !in.is_open() ){
339 fprintf(stderr,"Error: cannot open list file '%s'", filename);
340 return;
341 }
342
343 AliTPCCalPad *calPad=0x0;
344
345 TString sFile;
346 sFile.ReadFile(in);
347 in.close();
348
349 TObjArray *arrFileLine = sFile.Tokenize("\n");
350
351 TIter nextLine(arrFileLine);
352
353 TObjString *sObjLine=0x0;
354 while ( sObjLine = (TObjString*)nextLine() ){
355 TString sLine(sObjLine->GetString());
356
357 TObjArray *arrNextCol = sLine.Tokenize("\t");
358
359 TObjString *sObjType = (TObjString*)(arrNextCol->At(0));
360 TObjString *sObjFileName = (TObjString*)(arrNextCol->At(1));
361
362 if ( !sObjType || ! sObjFileName ) continue;
363 TString sType(sObjType->GetString());
364 TString sFileName(sObjFileName->GetString());
365 printf("%s\t%s\n",sType.Data(),sFileName.Data());
366
367 TFile *fIn = TFile::Open(sFileName);
368 if ( !fIn ){
369 fprintf(stderr,"File not found: '%s'", sFileName.Data());
370 continue;
371 }
372
373 if ( sType == "CE" ){
374 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
375
376 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
377 calPad->SetNameTitle("CETmean","CETmean");
378 calibObjects->Add(calPad);
379
380 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
381 calPad->SetNameTitle("CEQmean","CEQmean");
382 calibObjects->Add(calPad);
383
384 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
385 calPad->SetNameTitle("CETrms","CETrms");
386 calibObjects->Add(calPad);
387
388 } else if ( sType == "Pulser") {
0fe7645c 389 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
86df2b3a 390
391 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
392 calPad->SetNameTitle("PulserTmean","PulserTmean");
393 calibObjects->Add(calPad);
394
395 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
396 calPad->SetNameTitle("PulserQmean","PulserQmean");
397 calibObjects->Add(calPad);
398
399 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
400 calPad->SetNameTitle("PulserTrms","PulserTrms");
401 calibObjects->Add(calPad);
402
403 } else if ( sType == "Pedestals") {
404 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
405
406 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
407 calPad->SetNameTitle("Pedestals","Pedestals");
408 calibObjects->Add(calPad);
409
410 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
411 calPad->SetNameTitle("Noise","Noise");
412 calibObjects->Add(calPad);
413
414 } else {
415 fprintf(stderr,"Undefined Type: '%s'",sType.Data());
416
417 }
418 delete fIn;
419 }
420}
421
422
423
424void AliTPCcalibDB::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
425 //
426 // Write a tree with all available information
418bbcaf 427 // if mapFileName is specified, the Map information are also written to the tree
86df2b3a 428 // pads specified in outlierPad are not used for calculating statistics
429 // - the same function as AliTPCCalPad::MakeTree -
430 //
431 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
432
433 TObjArray* mapIROCs = 0;
434 TObjArray* mapOROCs = 0;
435 TVectorF *mapIROCArray = 0;
436 TVectorF *mapOROCArray = 0;
437 Int_t mapEntries = 0;
438 TString* mapNames = 0;
439
440 if (mapFileName) {
441 TFile mapFile(mapFileName, "read");
442
443 TList* listOfROCs = mapFile.GetListOfKeys();
444 mapEntries = listOfROCs->GetEntries()/2;
445 mapIROCs = new TObjArray(mapEntries*2);
446 mapOROCs = new TObjArray(mapEntries*2);
447 mapIROCArray = new TVectorF[mapEntries];
448 mapOROCArray = new TVectorF[mapEntries];
449
450 mapNames = new TString[mapEntries];
451 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
418bbcaf 452 TString nameROC(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
453 nameROC.Remove(nameROC.Length()-4, 4);
454 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "IROC").Data()), ivalue);
455 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((nameROC + "OROC").Data()), ivalue);
456 mapNames[ivalue].Append(nameROC);
86df2b3a 457 }
458
459 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
460 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
461 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
462
463 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
464 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
465 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
466 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
467 }
468
469 } // if (mapFileName)
470
471 TTreeSRedirector cstream(fileName);
472 Int_t arrayEntries = array->GetEntries();
473
474 TString* names = new TString[arrayEntries];
475 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
476 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
477
478 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
479 //
480 // get statistic for given sector
481 //
482 TVectorF median(arrayEntries);
483 TVectorF mean(arrayEntries);
484 TVectorF rms(arrayEntries);
485 TVectorF ltm(arrayEntries);
486 TVectorF ltmrms(arrayEntries);
487 TVectorF medianWithOut(arrayEntries);
488 TVectorF meanWithOut(arrayEntries);
489 TVectorF rmsWithOut(arrayEntries);
490 TVectorF ltmWithOut(arrayEntries);
491 TVectorF ltmrmsWithOut(arrayEntries);
492
493 TVectorF *vectorArray = new TVectorF[arrayEntries];
494 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
495 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
496
497 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
498 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
499 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
500 AliTPCCalROC* outlierROC = 0;
501 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
502 if (calROC) {
503 median[ivalue] = calROC->GetMedian();
504 mean[ivalue] = calROC->GetMean();
505 rms[ivalue] = calROC->GetRMS();
506 Double_t ltmrmsValue = 0;
507 ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
508 ltmrms[ivalue] = ltmrmsValue;
509 if (outlierROC) {
510 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
511 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
512 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
513 ltmrmsValue = 0;
514 ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
515 ltmrmsWithOut[ivalue] = ltmrmsValue;
516 }
517 }
518 else {
519 median[ivalue] = 0.;
520 mean[ivalue] = 0.;
521 rms[ivalue] = 0.;
522 ltm[ivalue] = 0.;
523 ltmrms[ivalue] = 0.;
524 medianWithOut[ivalue] = 0.;
525 meanWithOut[ivalue] = 0.;
526 rmsWithOut[ivalue] = 0.;
527 ltmWithOut[ivalue] = 0.;
528 ltmrmsWithOut[ivalue] = 0.;
529 }
530 }
531
532 //
533 // fill vectors of variable per pad
534 //
535 TVectorF *posArray = new TVectorF[8];
536 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
537 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
538
539 Float_t posG[3] = {0};
540 Float_t posL[3] = {0};
541 Int_t ichannel = 0;
542 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
543 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
544 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
545 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
546 posArray[0][ichannel] = irow;
547 posArray[1][ichannel] = ipad;
548 posArray[2][ichannel] = posL[0];
549 posArray[3][ichannel] = posL[1];
550 posArray[4][ichannel] = posG[0];
551 posArray[5][ichannel] = posG[1];
552 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
553 posArray[7][ichannel] = ichannel;
554
555 // loop over array containing AliTPCCalPads
556 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
557 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
558 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
559 if (calROC)
560 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
561 else
562 (vectorArray[ivalue])[ichannel] = 0;
563 }
564 ichannel++;
565 }
566 }
567
568 cstream << "calPads" <<
569 "sector=" << isector;
570
571 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
572 cstream << "calPads" <<
573 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
574 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
575 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
576 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
577 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
578 if (outlierPad) {
579 cstream << "calPads" <<
580 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
581 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
582 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
583 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
584 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
585 }
586 }
587
588 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
589 cstream << "calPads" <<
590 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
591 }
592
593 if (mapFileName) {
594 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
595 if (isector < 36)
596 cstream << "calPads" <<
597 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
598 else
599 cstream << "calPads" <<
600 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
601 }
602 }
603
604 cstream << "calPads" <<
605 "row.=" << &posArray[0] <<
606 "pad.=" << &posArray[1] <<
607 "lx.=" << &posArray[2] <<
608 "ly.=" << &posArray[3] <<
609 "gx.=" << &posArray[4] <<
610 "gy.=" << &posArray[5] <<
611 "rpad.=" << &posArray[6] <<
612 "channel.=" << &posArray[7];
613
614 cstream << "calPads" <<
615 "\n";
616
617 delete[] posArray;
618 delete[] vectorArray;
619 }
620
621
622 delete[] names;
623 if (mapFileName) {
624 delete mapIROCs;
625 delete mapOROCs;
626 delete[] mapIROCArray;
627 delete[] mapOROCArray;
628 delete[] mapNames;
629 }
630}
631