BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EmcBhaCalib.cxx
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2// File and Version Information:
3//
4// Description:
5// Class EmcBhaCalib - performs calibration of EMC Xtals with Bhabha
6// events and a Chi^2 fit, the resulting system of linear equations
7// of the fit is solved with the SLAP Algebra package
8//
9// Environment:
10// Software developed for the BESIII Detector at the BEPCII.
11//
12// Author List:
13// Chunxiu Liu
14//
15// Copyright Information:
16// Copyright (C) 2005 IHEP
17//
18//------------------------------------------------------------------------
19
20//-----------------------
21// This Class's Header --
22//-----------------------
24
25//-------------
26// C Headers --
27//-------------
28extern "C" {
29#include "slap/dlap.h"
30}
31#include <iostream>
32#include <fstream>
33#include <cmath>
34#include <cassert>
35#include <cstdlib>
36//-------------------------------
37// Collaborating Class Headers --
38//-------------------------------
39
40#include "GaudiKernel/MsgStream.h"
41#include "GaudiKernel/AlgFactory.h"
42#include "GaudiKernel/ISvcLocator.h"
43#include "GaudiKernel/SmartDataPtr.h"
44#include "GaudiKernel/IDataProviderSvc.h"
45#include "GaudiKernel/PropertyMgr.h"
46#include "GaudiKernel/DataObject.h"
47
48#include "CLHEP/Vector/ThreeVector.h"
49using namespace std;
50
51using CLHEP::Hep3Vector;
52
53//--------------------
54#include "GaudiKernel/MsgStream.h"
55
56
57#include "TROOT.h"
58#include "TFile.h"
59#include "TTree.h"
60#include "TH1F.h"
61#include "TObjArray.h"
62
63DECLARE_COMPONENT(EmcBhaCalib)
64
65DECLARE_COMPONENT(EmcBhaCalib)
66//----------------
67// Constructors --
68//----------------
69EmcBhaCalib::EmcBhaCalib(const std::string& name, ISvcLocator* pSvcLocator)
70 :Algorithm(name, pSvcLocator),
71 m_dirHitsCut(200),
72 m_convCrit(0.000001),
73 m_askForMatrixInversion(true),
74 m_fitPolynom(true), //no used now
75 m_writeToFile(false),
76 m_readDataFromDB(false),
77 m_equationSolver("GMR"),
78 m_fileExt(""),
79 m_fileDir("/home/besdata/public/liucx/Calib/"),
80 m_nrNonZeros(0),
81 m_nrXtalsEnoughHits(0),
82 m_runNumberFile("runnumbers.dat"),
83 m_MsgFlag(0)
84
85{
86
87 // Declare the properties
88 declareProperty("equationSolver", m_equationSolver);
89 declareProperty("dirHitsCut", m_dirHitsCut);
90 declareProperty("writeToFile", m_writeToFile);
91 declareProperty("fileExt", m_fileExt);
92 declareProperty("fileDir", m_fileDir);
93 declareProperty("readDataFromDB", m_readDataFromDB);
94 declareProperty("askForMatrixInversion", m_askForMatrixInversion);
95 declareProperty("fitPolynom", m_fitPolynom);
96 declareProperty("convCrit", m_convCrit);
97 declareProperty("runNumberFile", m_runNumberFile);
98 declareProperty("MsgFlag", m_MsgFlag);
99
100 m_calibConst = new double[6240];
101 m_calibConstUnred = new double[6240];
102 m_absoluteConstants = new double[6240];
103 m_oldConstants = new double[6240];
104
105
106 //ntuple
107 m_tuple1=0;
108
109
110
111}
112
113//--------------
114// Destructor --
115//--------------
117 if ( 0 != m_calibConst) {
118 delete [] m_calibConst;
119 m_calibConst = 0;
120 }
121 if ( 0 != m_calibConstUnred) {
122 delete [] m_calibConstUnred;
123 m_calibConstUnred = 0;
124 }
125 if ( 0 != m_absoluteConstants) {
126 delete [] m_absoluteConstants;
127 m_absoluteConstants = 0;
128 }
129 if ( 0 != m_oldConstants) {
130 delete [] m_oldConstants;
131 m_oldConstants = 0;
132 }
133 if ( 0 != m_myCalibData) {
134 delete m_myCalibData;
135 m_myCalibData = 0;
136 }
137
138}
139
140// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
142
143 MsgStream log(msgSvc(), name());
144 log << MSG::INFO << "in initialize()" << endreq;
145
146 m_myCalibData = new EmcBhaCalibData(6240, m_MsgFlag);
147
148 m_calibrationSuccessful = false;
149
150 StatusCode status1;
151
152 NTuplePtr nt1(ntupleSvc(),"FILE302/n1");
153 if ( nt1 ) m_tuple1 = nt1;
154 else {
155 m_tuple1=ntupleSvc()->book("FILE302/n1",CLID_ColumnWiseTuple,"Calib");
156 if( m_tuple1 ) {
157
158 status1 = m_tuple1->addItem ("ixtal",ixtal);
159 status1 = m_tuple1->addItem ("gi0",gi0);
160 status1 = m_tuple1->addItem ("calibConst",calibConst);
161 status1 = m_tuple1->addItem ("err",err);
162 status1 = m_tuple1->addItem ("nhitxtal",nhitxtal);
163
164 }
165 else { // did not manage to book the N tuple....
166 log << MSG::ERROR <<"Cannot book N-tuple:" << long(m_tuple1) << endmsg;
167 return StatusCode::FAILURE;
168 }
169 }
170
171 // use EmcCalibConstSvc
172 StatusCode scCalib;
173 scCalib = Gaudi::svcLocator() -> service("EmcCalibConstSvc", m_emcCalibConstSvc);
174 if( scCalib != StatusCode::SUCCESS){
175 log << MSG::ERROR << "can not use EmcCalibConstSvc" << endreq;
176 }
177 else {
178 std::cout << "Test EmcCalibConstSvc DigiCalibConst(0)= "
179 << m_emcCalibConstSvc -> getDigiCalibConst(0) << std::endl;
180 }
181
182
183 //std::cout<<6239<<"\t"<< m_emcCalibConstSvc->getCrystalEmaxData(6239) <<endl;
184
185
186 //init starting values for calibration constants from file or set to 1
187 initCalibConst();
188
189 //digiConstCor();
190
191 return StatusCode::SUCCESS;
192}
193
194// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
196
197 MsgStream log(msgSvc(), name());
198 log << MSG::DEBUG << "in execute()" << endreq;
199 //read in accumulated matrix/matrices
200 if ( m_readDataFromDB ) {
201 m_calibrationSuccessful = readInFromDB();
202
203 log << MSG::INFO << "read in accumulated matrix from DB"<<endreq;
204
205 } else {
206 m_calibrationSuccessful = readInFromFile();
207
208 log << MSG::INFO <<"read in accumulated matrix from file"<<endreq;
209
210 }
211
212 if ( m_calibrationSuccessful == true ) {
213
214 //ask first if to do a matrix inversion
215 if (m_askForMatrixInversion==true) {
216
217 m_calibrationSuccessful = false;
218 log << MSG::INFO << " Well, we have "
219 << m_nrXtalsEnoughHits << " crystals whith more "
220 << "than " << m_dirHitsCut
221 << " direct hits. Do you want do "
222 << "a matrix inversion ? [N]" << endreq;
223
224
225 m_calibrationSuccessful = true;
226
227 }
228
229 if ( m_calibrationSuccessful == true ) {
230
231 //write added matrix and vector to files
232 if ( m_writeToFile == true) {
233 writeOut();
234 }
235 //cout <<"writeout"<<endl;
236 m_myCalibData->printVec(2);
237
238 //reduce to continious array of only non zeros elements for SLAP
239 if ( m_calibrationSuccessful = m_myCalibData->reduce() ) {
240
241 cout<<"m_calibrationSuccessful="<<m_calibrationSuccessful<<endl;
242
243 if ( m_myCalibData->nXtalsHit() != m_myCalibData->nXtals() ){
244 log << MSG::INFO << " Reduced number of Xtals for calibration: "
245 << m_myCalibData->nXtalsHit()
246 << endreq;
247 }
248 cout<<"m_myCalibData->nXtalsHit()="<<m_myCalibData->nXtalsHit()
249 <<"m_myCalibData->nXtals()="<<m_myCalibData->nXtals()<<endl;
250 m_myCalibData->printVec(10);
251
252 //solve matrix equation
253 if ( m_calibrationSuccessful = solveEqua() ) {
254
255 for (int i=0; i<m_myCalibData->nXtalsHit(); i++){
256 //fill an array of constants for storage in database
257 m_calibConstUnred[m_myCalibData->xtalInd(i)]
258 = m_calibConst[i];
259
260 // if (m_myCalibData->xtalHitsDir(i)>10)
261 // cout<<"Index,calibconst="<<" "<<i <<" "<<m_myCalibData->xtalInd(i)
262 // <<" "<<m_calibConstUnred[m_myCalibData->xtalInd(i)]
263 // <<" "<<m_calibConstUnred[m_myCalibData->xtalInd(i)]*
264 // m_oldConstants[m_myCalibData->xtalInd(i)]<<endl;
265
266 }
267 //do validation, fit polynom if necessary, fill CalList
268 prepareConstants();
269
270 //if wanted write constants to plain ASCII file
271 if ( m_writeToFile==true){
272 writeOutConst();
273 }
274
275 ntupleOut();
276
277 }
278
279 } else {
280 log << MSG::WARNING << " Reduction of the Emc Bhabha calibration matrix FAILED !"
281 << endl
282 << " Will NOT perform Emc Bhabha calibration !"
283 << endreq;
284 return( StatusCode::FAILURE);
285 }
286 }
287 } else {
288 log << MSG::WARNING << " ERROR in reading in Emc Bhabha calibration data !"
289 << endl
290 << " Will NOT perform Emc Bhabha calibration !"
291 << endreq;
292 return( StatusCode::FAILURE);
293 }
294
295 return StatusCode::SUCCESS;
296}
297
298// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
300
301 MsgStream log(msgSvc(), name());
302
303
304 log << MSG::INFO << "in endRun()" << endreq;
305
306
307 return StatusCode::SUCCESS;
308}
309
310
311void
313
314 MsgStream log(msgSvc(), name());
315
316 log << MSG::INFO<< "Performs the Chi square fit of the accumulated "
317 << endreq;
318}
319
320//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
321void
322EmcBhaCalib::initCalibConst( ) {
323
324
325 MsgStream log(msgSvc(), name());
326
327 for ( int i = 0; i< m_myCalibData->nXtals(); i++ ) {
328 m_calibConst[i] = 1.;
329 }
330
331 ifstream inConst("EmcCalib.Const");
332 if ( !inConst ) {
333 log << MSG::VERBOSE
334 << " No starting values for calibration constants found ! "
335 << "(EmcCalib.Const)"
336 << endl
337 << "Set them to 1.0 ! " << endreq;
338 } else {
339 log << MSG::VERBOSE << " Read in starting values of calibration constants !"
340 << endreq;
341
342 int nConst;
343 inConst >> nConst;
344 log << MSG::VERBOSE << " Have starting values for "
345 << nConst << " Xtals !" << endl
346 << "Set the others to 1.0 ! " << endmsg;
347
348 int numberXtal;
349 for ( int i = 0; i< nConst; i++ ) {
350 inConst >> numberXtal >> m_calibConst[numberXtal];
351 }
352 }
353
354 int nConstEmc;
355
356 nConstEmc= m_emcCalibConstSvc -> getDigiCalibConstNo() ;
357
358 if ( nConstEmc!=6240) cout<<"number of calibconst="<< nConstEmc<<endl;
359
360 // log << MSG::VERBOSE << " Have starting values for "
361 // << nConstEmc << " Xtals !" << endl
362 // << "Set the others to 1.0 ! " << endmsg;
363
364 for ( int i = 0; i< nConstEmc; i++ ) {
365 //cout<<i<<" "
366 //<<m_emcCalibConstSvc->getThetaIndex(i)<<" "
367 //<<m_emcCalibConstSvc->getPhiIndex(i)
368 //<<" "<<m_emcCalibConstSvc->getEmcID(i)<<endl;
369 m_calibConst[i] = m_emcCalibConstSvc -> getDigiCalibConst(i);
370 //m_calibConst[i] =5.0;
371 m_oldConstants[i]=m_emcCalibConstSvc -> getDigiCalibConst(i);
372
373 // initialize m_absoluteConstants as m_oldConstants.
374 //m_absoluteConstants[i] =m_emcCalibConstSvc -> getDigiCalibConst(i);
375 m_absoluteConstants[i] =1.0;
376 // initialize m_calibConstUnred as 1.
377 m_calibConstUnred[i] =1.0;
378 }
379
380
381}
382
383
384//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
385bool
386EmcBhaCalib::solveEqua() {
387
388 MsgStream log(msgSvc(), name());
389 //-----------------------------------------------------
390 // solve system of equations with SLAP package
391 //-----------------------------------------------------
392 log << MSG::VERBOSE << " Solve Matrix equation with method "
393 << m_equationSolver
394 << endl
395 << "Xtals hit: " << m_myCalibData->nXtalsHit()
396 << endl
397 << "Nr of non Zeros: " << m_nrNonZeros << endreq;
398
399 //input parameters for SLAP
400 long int isym = 1; //matrix M is symmetric
401 long int nsave = 20;
402 long int itol = 0; //type of convergence criterion
403 double convCriterion = m_convCrit; //convergence crtiterion
404 long int maxIter = 1000; //maximum number of iterations
405 long int errUnit = 6; //unit on which to write error estimation
406 // at each iteration
407
408 //working arrays needed by SLAP
409 long int lengthDoubleWork;
410 double* doubleWork;
411 long int lengthIntWork;
412 long int* intWork;
413
414 //output parameters of slap
415 long int iters=0; //number of iterations required
416 double errorv=1000; //error estimate
417 long int errorFlag=9999;
418
419 //sparse Ax=b solver.
420 //uses the generalized minimum residual method
421 //The preconditioner is diagonal scaling.
422 if ( m_equationSolver == "GMR" ) {
423
424 cout<<m_equationSolver<<endl;
425
426 cout<<"errorFlag="<<errorFlag<<endl;
427
428 //workspaces
429 lengthDoubleWork = (1 + m_myCalibData->nXtals()*(nsave+7)
430 + nsave*(nsave+3))
431 + 50000;
432 doubleWork = new double[lengthDoubleWork];
433 lengthIntWork = 50;
434 intWork = new long int [lengthIntWork];
435
436 dsdgmr_ ( &(m_myCalibData->nXtalsHit()),
437 m_myCalibData->getVectorR(),
438 m_calibConst,
439 &m_nrNonZeros,
440 m_myCalibData->getMatrixM()->row(),
441 m_myCalibData->getMatrixM()->column(),
442 m_myCalibData->getMatrixM()->matrix(),
443 &isym,
444 &nsave,
445 &itol,
446 &convCriterion,
447 &maxIter,
448 &iters,
449 &errorv,
450 &errorFlag,
451 &errUnit,
452 doubleWork,
453 &lengthDoubleWork,
454 intWork,
455 &lengthIntWork
456 );
457 cout<<"errorFlag="<<errorFlag<<endl;
458
459
460 log << MSG::VERBOSE << " Error flag: " << errorFlag << endreq;
461 if ( errorFlag < 0 ) errorFlag = labs(errorFlag) + 2;
462 switch ( errorFlag ) {
463 case 0: log << MSG::VERBOSE << " all went well"
464 << endreq; break;
465 case 1: log << MSG::ERROR << " insufficient storage allocated for _doubleWork or _intWork"
466 << endreq; break;
467 case 2: log << MSG::ERROR << " method failed to reduce norm of current residual"
468 << endreq; break;
469 case 3: log << MSG::ERROR << " insufficient length for _doubleWork"
470 << endreq; break;
471 case 4: log << MSG::ERROR << " inconsistent _itol and values"
472 << endreq; break;
473 default: log << MSG::ERROR << " unknown flag" << endreq;
474 }
475 log << MSG::VERBOSE << " Integer workspace used = " << intWork[8] << endl
476 << " Double workspace used = " << intWork[9] << endreq;
477 }
478
479 //Routine to solve a symmetric positive definite linear
480 //system Ax = b using the Preconditioned Conjugate
481 //Gradient method. The preconditioner is diagonal scaling.
482 if ( m_equationSolver == "PCG" ) {
483
484 cout<<m_equationSolver<<endl;
485
486 itol = 1;
487
488 //workspaces
489 lengthDoubleWork = 5 *m_myCalibData->nXtals() + 50000;
490 doubleWork = new double[lengthDoubleWork];
491 lengthIntWork = 50;
492 intWork = new long int [lengthIntWork];
493
494 dsdcg_( &(m_myCalibData->nXtalsHit()),
495 m_myCalibData->getVectorR(),
496 m_calibConst,
497 &m_nrNonZeros,
498 m_myCalibData->getMatrixM()->row(),
499 m_myCalibData->getMatrixM()->column(),
500 m_myCalibData->getMatrixM()->matrix(),
501 &isym,
502 &itol,
503 &convCriterion,
504 &maxIter,
505 &iters,
506 &errorv,
507 &errorFlag,
508 &errUnit,
509 doubleWork,
510 &lengthDoubleWork,
511 intWork,
512 &lengthIntWork
513 );
514
515 switch ( errorFlag ) {
516 case 0: log << MSG::VERBOSE << "all went well" << endreq; break;
517 case 1: log << MSG::ERROR << " insufficient storage allocated for WORK or IWORK"
518 << endreq; break;
519 case 2: log << MSG::ERROR << " method failed to to converge in maxIter steps."
520 << endreq; break;
521 case 3:log << MSG::ERROR << " Error in user input. Check input value of _nXtal,_itol."
522 << endreq; break;
523 case 4:log << MSG::ERROR << " User error tolerance set too tight. "
524 << "Reset to 500.0*D1MACH(3). Iteration proceeded."
525 << endreq; break;
526 case 5:log << MSG::ERROR << " Preconditioning matrix, M, is not Positive Definite. "
527 << endreq; break;
528 case 6: log << MSG::ERROR << " Matrix A is not Positive Definite."
529 << endreq; break;
530 default: log << MSG::ERROR << " unknown flag" << endreq;
531 }
532 log << MSG::VERBOSE << " Integer workspace used = " << intWork[9] << endl
533 << "Double workspace used = " << intWork[10] << endreq;
534 }
535
536 log << MSG::VERBOSE << "------ Calibration fit statistics ------- " << endl
537 << "maximum number of iterations =" << maxIter << endl
538 << " number of iterations =" << iters << endl
539 << "error estimate of error in final solution ="
540 << errorv << endreq;
541
542 if ( 0 != doubleWork) delete [] doubleWork;
543 if ( 0 != intWork) delete [] intWork;
544
545 if ( errorFlag != 0 ) return false;
546 else return true;
547
548 return true;
549}
550
551//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
552void
553EmcBhaCalib::writeOut() {
554
555 ofstream vectorOut;
556 std::string vectorFile = m_fileDir;
557 vectorFile += m_fileExt;
558 vectorFile += "allCalibVector.dat";
559 vectorOut.open(vectorFile.c_str());
560
561 ofstream matrixOut;
562 std::string matrixFile = m_fileDir;
563 matrixFile += m_fileExt;
564 matrixFile += "allCalibMatrix.dat";
565 matrixOut.open(matrixFile.c_str());
566
567 MsgStream log(msgSvc(), name());
568
569 log << MSG::VERBOSE << " Write out files "
570 << vectorFile << " "
571 << matrixFile
572 << endreq;
573
574 m_myCalibData->writeOut(matrixOut,vectorOut);
575
576 vectorOut.close();
577 matrixOut.close();
578
579}
580
581//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
582void
583EmcBhaCalib::writeOutConst() {
584
585 ofstream constOut;
586
587 std::string constFile = m_fileDir;
588 constFile += m_fileExt;
589 constFile += "CalibConst.dat";
590
591 constOut.open(constFile.c_str());
592
593 constOut << "#crystalIndex relative-constant absolute-constant" << endl;
594
595 int chanIndex;
596
597 //output constants to file
598 for ( int i=0; i < m_myCalibData->nXtalsHit(); i++) {
599
600 chanIndex=m_myCalibData->xtalInd(i); //xtalInd(i) array of xtal indices
601
602 //---- get the new absolute constant ----
603 //set it to 0 if we have not enough hits -> we'll keep the old constant
604 if ( m_myCalibData->xtalHitsDir(i) < m_dirHitsCut )
605 m_absoluteConstants[chanIndex] = m_oldConstants[chanIndex];
606 else
607 m_absoluteConstants[chanIndex] =
608 m_oldConstants[chanIndex] * m_calibConstUnred[chanIndex];
609
610 }
611
612 //output constants to file
613 for ( int i=0; i < m_myCalibData->nXtals(); i++) {
614
615 long Xtal_Index = i;
616 if(m_calibConstUnred[Xtal_Index]>0){
617 constOut << Xtal_Index << " "
618 << m_calibConstUnred[Xtal_Index] << " "
619 << m_oldConstants[Xtal_Index] << " "
620 << m_absoluteConstants[Xtal_Index]
621 << endl;
622
623
624 }
625 }
626 constOut.close();
627
628}
629
630//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
631bool
632EmcBhaCalib::readInFromFile() {
633
634 MsgStream log(msgSvc(), name());
635
636 log << MSG::INFO << " Read Bhabha calibration data from files !"
637 << endreq;
638
639 std::string runNumberFile = m_fileDir;
640 runNumberFile += m_fileExt;
641 runNumberFile += m_runNumberFile;
642
643 bool success = false;
644 log << MSG::INFO << " Get file names from input file : "
645 << runNumberFile
646 << endreq;
647
648 std::string vectorFile;
649 std::string matrixFile;
650
651 log << MSG::INFO << "Runnumber non-zeros xtals"
652 << endreq;
653
654 ifstream datafile(runNumberFile.c_str());
655
656 //cout<<"datafile="<<runNumberFile.c_str()<<""<<datafile.good()<<endl;
657
658 if (datafile.good() > 0 ) {
659
660 while( !datafile.eof() ) {
661
662 ifstream vectorIn;
663 ifstream matrixIn;
664
665 std::string num;
666 datafile >> num;
667
668 if ( num.length() > 0 ) {
669
670 vectorFile = m_fileDir;
671 vectorFile += m_fileExt;
672 vectorFile += num;
673 vectorFile += "CalibVector.dat";
674
675 matrixFile = m_fileDir;
676 matrixFile += m_fileExt;
677 matrixFile += num;
678 matrixFile += "CalibMatrix.dat";
679
680 vectorIn.open(vectorFile.c_str());
681 matrixIn.open(matrixFile.c_str());
682
683 if ( vectorIn.good() > 0 && matrixIn.good() > 0 ) {
684
685 m_myCalibData->readIn(matrixIn,vectorIn);
686
687 log << MSG::INFO << num
688 << " "
689 << m_myCalibData->getMatrixM()->num_nonZeros()
690 << " "
691 << m_myCalibData->nXtalsHit()
692 << endreq;
693
694 success = true;
695
696 }
697 else {
698 if ( !vectorIn )
699 log << MSG::ERROR << " ERROR: Vector input file "
700 << vectorFile
701 << " doesn't exist !" << endreq;
702 if ( !matrixIn )
703 log << MSG::ERROR << " ERROR: Matrix input file "
704 << matrixFile
705 << " doesn't exist !" << endreq;
706 }
707 vectorIn.close();
708 matrixIn.close();
709 }
710 }
711 }
712
713 if ( success == true) {
714
715 for (int i=0; i < m_myCalibData->nXtals(); i++) {
716
717 if ( m_myCalibData->xtalHitsDir(i) > m_dirHitsCut )
718 {
719 m_nrXtalsEnoughHits++;
720 }
721 }
722 m_nrNonZeros = m_myCalibData->getMatrixM()->num_nonZeros();
723 log << MSG::INFO<< " Final matrix : "
724 << "Number of non zero elements: " << m_nrNonZeros
725 << " "
726 << m_myCalibData->nXtalsHit()
727 << endreq;
728
729 }
730
731
732 return success;
733}
734
735//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
736bool
737EmcBhaCalib::readInFromDB() {
738
739 bool success = true;
740
741 return success;
742}
743
744//* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
745bool
746EmcBhaCalib::prepareConstants() {
747
748 bool successfull=false;
749
750 //the old constant
751 //float theCalConst = 1.0;
752 int chanIndex;
753
754 for ( int i = 0; i< m_myCalibData->nXtalsHit(); i++ ) {
755
756 chanIndex=m_myCalibData->xtalInd(i); //xtalInd(i) array of xtal indices
757
758 //////////deal with the crystal of ixtal=689, because of ADC is very small;
759 // m_calibConstUnred[chanIndex]=1;
760 //if (chanIndex==689) m_calibConstUnred[chanIndex]=1.47;
761 //if (chanIndex==689)
762 // m_oldConstants[chanIndex]=m_oldConstants[chanIndex]*1.47;
763 //---- get the new absolute constant ----
764 //set it to 0 if we have not enough hits -> we'll keep the old constant
765 if ( m_myCalibData->xtalHitsDir(i) < m_dirHitsCut )
766 m_absoluteConstants[chanIndex] = m_oldConstants[chanIndex];
767 else
768 m_absoluteConstants[chanIndex] =
769 m_oldConstants[chanIndex] * m_calibConstUnred[chanIndex];
770 }
771
772
773 double DigiConst[6240];
774 int IxtalNumber[6240];
775
776 for(int ind=0; ind<m_myCalibData->nXtals(); ind++ ) {
777
778 DigiConst[ind]=m_absoluteConstants[ind];
779 // cout<<"ind="<<ind<<"\t"<< DigiConst[ind]<<endl;
780 }
781
782 int ixtal1,ixtal2,ixtal3;
783 ixtal1=m_emcCalibConstSvc->getIndex(1,20,26);
784 ixtal2=m_emcCalibConstSvc->getIndex(1,23,26);
785 ixtal3=m_emcCalibConstSvc->getIndex(1,24,26);
786 //cout<<ixtal1<<" "<<ixtal2<<" "<<ixtal3<<" "<<endl;
787 for(int ind=0; ind<m_myCalibData->nXtals(); ind++ ) {
788
789 IxtalNumber[ind]=-1;
790 /*
791 if (ind==ixtal1) IxtalNumber[ind]=ixtal3;
792 if (ind==ixtal2) IxtalNumber[ind]=ixtal1;
793 if (ind==ixtal3) IxtalNumber[ind]=ixtal2;
794 */
795 }
796
797 TFile fconst("EmcCalibConst.root", "recreate");
798
799 //define tree fill the absolute digicalibconsts into the root-file
800 TTree* constr= new TTree("DigiCalibConst","DigiCalibConst");
801 constr->Branch("DigiCalibConst",DigiConst,"DigiConst[6240]/D");
802 constr->Branch("IxtalNumber",IxtalNumber,"IxtalNumber[6240]/I");
803
804 constr->Fill();
805 /*
806
807 double aa[6240];
808 int Iaa[6240];
809 constr->SetBranchAddress("DigiCalibConst", aa);
810 constr->SetBranchAddress("IxtalNumber", Iaa);
811 constr->GetEntry(0);
812 cout<<Iaa[10]<<endl;
813 cout<<aa[10]<<endl;
814
815 cout<<constr->GetEntry(0)<<endl;
816
817 */
818
819 constr->Write();
820
821 delete constr;
822
823 fconst.Close();
824
825 // end tree
826
827 successfull=true;
828 return successfull;
829
830}
831
832void
833EmcBhaCalib::ntupleOut(){
834
835
836 for (int i=0; i < m_myCalibData->nXtalsHit(); i++) {
837 int xtalIndex=m_myCalibData->xtalInd(i);
838
839 gi0 = m_oldConstants[xtalIndex];
840
841 if (gi0<0) cout<<"gi0="<<gi0<<endl;
842 err = (m_calibConstUnred[xtalIndex]-gi0)/gi0;
843 calibConst=m_calibConstUnred[xtalIndex];
844 ixtal=xtalIndex;
845 nhitxtal=m_myCalibData->xtalHitsDir(i); ///
846 m_tuple1->write();
847 }
848
849}
850
851
852void
853EmcBhaCalib::printInfo(std::string fileName) {
854
855 ofstream xtalHitsDirOut;
856 std::string xtalHitsDirFile = m_fileDir;
857 xtalHitsDirFile += m_fileExt;
858 xtalHitsDirFile += fileName;
859 xtalHitsDirOut.open(xtalHitsDirFile.c_str());
860
861 //nXtalsHit() is number of xtals hit
862 for (int i=0; i < m_myCalibData->nXtalsHit(); i++) {
863 if ( m_myCalibData->xtalHitsDir(i) > m_dirHitsCut )
864 {
865 int chanindex=m_myCalibData->xtalInd(i);
866 xtalHitsDirOut<<chanindex<<endl;
867 }
868 }
869 xtalHitsDirOut.close();
870
871}
872
873void
874EmcBhaCalib::digiConstCor(){
875
876 // ofstream Fuout;
877 //std::string Fuoutfile;
878
879 //Fuoutfile="emccalibconst.txt";
880 //Fuout.open(Fuoutfile.c_str());
881
882 double digiConst[6240];
883
884 for(int ind=0; ind<m_myCalibData->nXtals(); ind++ ) {
885
886 digiConst[ind]=m_oldConstants[ind];
887 // Fuout<<m_oldConstants[ind]<<endl;
888 }
889
890 // Fuout.close();
891
892 ifstream ExpEneIn;
893 std::string ExpEneFile;
894 ExpEneFile = "cor.conf";
895 ExpEneIn.open(ExpEneFile.c_str());
896
897 double ExpEne[56];
898 double Exp_barrel[44],Exp_east[6],Exp_west[6];
899
900 for(int i=0;i<56;i++) {
901
902 ExpEneIn>>ExpEne[i];
903 ExpEne[i]=ExpEne[i]*1.843;
904 if (i<6) Exp_east[i]=ExpEne[i];
905 if (i>=6&&i<50) Exp_barrel[i-6]=ExpEne[i];
906 if (i>=50) Exp_west[55-i]=ExpEne[i];
907
908 //cout<<ExpEne[i]<<endl;
909 }
910 ExpEneIn.close();
911
912 ifstream EDepEneIn;
913 std::string EDepEneFile = "allxtal.dat";
914 EDepEneIn.open(EDepEneFile.c_str());
915
916 double CorDigiConst[6240];
917 int ipart,ithe,iphi;
918 double peak,sigma;
919 double coeff;
920
921 for(int i=0;i<6240;i++) {
922
923 EDepEneIn>>ipart>>ithe>>iphi>>peak>>sigma;
924 int ix=m_emcCalibConstSvc->getIndex(ipart,ithe,iphi);
925
926 if (ipart==0) coeff=peak/Exp_east[ithe];
927 if (ipart==1) coeff=peak/Exp_barrel[ithe];
928 if (ipart==2) coeff=peak/Exp_west[ithe];
929 cout<<coeff<<endl;
930 //CorDigiConst[ix]=digiConst[ix]/coeff;
931 CorDigiConst[ix]=coeff;
932 }
933 EDepEneIn.close();
934
935
936 //define tree fill the absolute digicalibconsts into the root-file
937
938 TTree* constr= new TTree("DigiCalibConst","DigiCalibConst");
939 constr->Branch("DigiCalibConst",CorDigiConst,"CorDigiConst[6240]/D");
940 constr->Fill();
941
942 double aa[6240];
943 constr->SetBranchAddress("DigiCalibConst", aa);
944 constr->GetEntry(0);
945
946 TFile fconst("EmcCalibConstCor.root", "recreate");
947 constr->Write();
948
949 fconst.Close();
950 delete constr;
951
952
953}
TTree * sigma
int dsdcg_(const int *n, const double *b, double *x, const long *nelt, int *ia, int *ja, double *a, const long *isym, const long *itol, const double *tol, const long *itmax, long *iter, double *err, long *ierr, const long *iunit, double *rwork, const long *lenw, long *iwork, const long *leniw)
int dsdgmr_(const int *n, const double *b, double *x, const long *nelt, int *ia, int *ja, double *a, const long *isym, long *nsave, const long *itol, const double *tol, const long *itmax, long *iter, double *err, long *ierr, const long *iunit, double *rwork, const long *lenw, long *iwork, const long *leniw)
NTuple::Tuple * m_tuple1
Definition: MdcHistItem.h:45
INTupleSvc * ntupleSvc()
IMessageSvc * msgSvc()
int & xtalHitsDir(int ind)
double * getVectorR()
void printVec(int number)
int xtalInd(int ind)
EmcLSSMatrix * getMatrixM()
void writeOut(ostream &OutM, ostream &outV)
void readIn(istream &InM, istream &InV)
StatusCode execute()
virtual void help()
StatusCode finalize()
virtual ~EmcBhaCalib()
StatusCode initialize()
int * row() const
Definition: EmcLSSMatrix.h:84
long int num_nonZeros()
int * column(const int &rowind=0) const
double * matrix(const int &rowind=0) const
virtual int getIndex(unsigned int PartId, unsigned int ThetaIndex, unsigned int PhiIndex) const =0
#define ix(i)
std::ifstream ifstream
Definition: bpkt_streams.h:44
std::ofstream ofstream
Definition: bpkt_streams.h:42