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